change gmx_bool to bool in gmxpreprocess
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / readir.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef GMX_GMXPREPROCESS_READIR_H
39 #define GMX_GMXPREPROCESS_READIR_H
40
41 #include "gromacs/fileio/readinp.h"
42 #include "gromacs/gmxpreprocess/grompp-impl.h"
43
44 namespace gmx
45 {
46 class MDModules;
47 }
48
49 struct gmx_groups_t;
50 struct gmx_mtop_t;
51 struct gmx_output_env_t;
52 struct pull_params_t;
53 struct pull_t;
54 struct t_grpopts;
55 struct t_inpfile;
56 struct t_inputrec;
57 struct t_rot;
58 struct warninp;
59 typedef warninp *warninp_t;
60
61 enum {
62     eshNONE, eshHBONDS, eshALLBONDS, eshHANGLES, eshALLANGLES, eshNR
63 };
64
65 enum {
66     ecouplamVDWQ, ecouplamVDW, ecouplamQ, ecouplamNONE, ecouplamNR
67 };
68
69 struct t_gromppopts
70 {
71     int      warnings;
72     int      nshake;
73     char    *include;
74     char    *define;
75     bool     bGenVel;
76     bool     bGenPairs;
77     real     tempi;
78     int      seed;
79     bool     bOrire;
80     bool     bMorse;
81     char    *wall_atomtype[2];
82     char    *couple_moltype;
83     int      couple_lam0;
84     int      couple_lam1;
85     bool     bCoupleIntra;
86 };
87
88 /*! \brief Initialise object to hold strings parsed from an .mdp file */
89 void init_inputrec_strings();
90
91 /*! \brief Clean up object that holds strings parsed from an .mdp file */
92 void done_inputrec_strings();
93
94 void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
95               warninp_t wi);
96 /* Validate inputrec data.
97  * Fatal errors will be added to nerror.
98  */
99 int search_string(const char *s, int ng, char *gn[]);
100 /* Returns the index of string s in the index groups */
101
102 void double_check(t_inputrec *ir, matrix box,
103                   bool bHasNormalConstraints,
104                   bool bHasAnyConstraints,
105                   warninp_t wi);
106 /* Do more checks */
107
108 void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
109                   warninp_t wi);
110 /* Do even more checks */
111
112 void check_chargegroup_radii(const gmx_mtop_t *mtop, const t_inputrec *ir,
113                              rvec *x,
114                              warninp_t wi);
115 /* Even more checks, charge group radii vs. cut-off's only. */
116
117 void get_ir(const char *mdparin, const char *mdparout,
118             gmx::MDModules *mdModules, t_inputrec *ir, t_gromppopts *opts,
119             WriteMdpHeader writeMdpHeader, warninp_t wi);
120 /* Read the input file, and retrieve data for inputrec.
121  * More data are read, but the are only evaluated when the next
122  * function is called. Also prints the input file back to mdparout.
123  */
124
125 void do_index(const char* mdparin,
126               const char *ndx,
127               gmx_mtop_t *mtop,
128               bool        bVerbose,
129               t_inputrec *ir,
130               warninp_t   wi);
131 /* Read the index file and assign grp numbers to atoms.
132  */
133
134 /* Routines In readpull.c */
135
136 char **read_pullparams(std::vector<t_inpfile> *inp,
137                        pull_params_t          *pull,
138                        warninp_t               wi);
139 /* Reads the pull parameters, returns a list of the pull group names */
140
141 void make_pull_groups(pull_params_t *pull,
142                       char **pgnames,
143                       const t_blocka *grps, char **gnames);
144 /* Process the pull group parameters after reading the index groups */
145
146 void make_pull_coords(pull_params_t *pull);
147 /* Process the pull coordinates after reading the pull groups */
148
149 pull_t *set_pull_init(t_inputrec *ir, const gmx_mtop_t *mtop,
150                       rvec *x, matrix box, real lambda);
151 /* Prints the initial pull group distances in x.
152  * If requested, adds the current distance to the initial reference location.
153  * Returns the pull_t pull work struct. This should be passed to finish_pull()
154  * after all modules have registered their external potentials, if present.
155  */
156
157 int str_nelem(const char *str, int maxptr, char *ptr[]);
158 /* helper function from readir.c to convert strings */
159
160 char **read_rotparams(std::vector<t_inpfile> *inp, t_rot *rot, warninp_t wi);
161 /* Reads enforced rotation parameters, returns a list of the rot group names */
162
163 void make_rotation_groups(t_rot *rot, char **rotgnames,
164                           t_blocka *grps, char **gnames);
165 /* Process the rotation parameters after reading the index groups */
166
167 void set_reference_positions(t_rot *rot, rvec *x, matrix box,
168                              const char *fn, bool bSet, warninp_t wi);
169
170 #endif