3c2a08577800ac7b7787565fdd18c97c6e875762
[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,2019, 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/math/vectypes.h"
43 #include "gromacs/utility/real.h"
44
45 namespace gmx
46 {
47 class MDModules;
48 }
49
50 struct gmx_mtop_t;
51 struct gmx_output_env_t;
52 struct pull_params_t;
53 struct pull_t;
54 struct t_blocka;
55 struct t_grpopts;
56 struct t_inpfile;
57 struct t_inputrec;
58 struct t_rot;
59 struct warninp;
60 typedef warninp *warninp_t;
61
62 enum {
63     eshNONE, eshHBONDS, eshALLBONDS, eshHANGLES, eshALLANGLES, eshNR
64 };
65
66 enum {
67     ecouplamVDWQ, ecouplamVDW, ecouplamQ, ecouplamNONE, ecouplamNR
68 };
69
70 struct t_gromppopts
71 {
72     int      warnings;
73     int      nshake;
74     char    *include;
75     char    *define;
76     bool     bGenVel;
77     bool     bGenPairs;
78     real     tempi;
79     int      seed;
80     bool     bOrire;
81     bool     bMorse;
82     char    *wall_atomtype[2];
83     char    *couple_moltype;
84     int      couple_lam0;
85     int      couple_lam1;
86     bool     bCoupleIntra;
87 };
88
89 /*! \brief Initialise object to hold strings parsed from an .mdp file */
90 void init_inputrec_strings();
91
92 /*! \brief Clean up object that holds strings parsed from an .mdp file */
93 void done_inputrec_strings();
94
95 void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
96               warninp_t wi);
97 /* Validate inputrec data.
98  * Fatal errors will be added to nerror.
99  */
100 int search_string(const char *s, int ng, char *gn[]);
101 /* Returns the index of string s in the index groups */
102
103 void double_check(t_inputrec *ir, matrix box,
104                   bool bHasNormalConstraints,
105                   bool bHasAnyConstraints,
106                   warninp_t wi);
107 /* Do more checks */
108
109 void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
110                   warninp_t wi);
111 /* Do even more checks */
112
113 void check_chargegroup_radii(const gmx_mtop_t *mtop, const t_inputrec *ir,
114                              rvec *x,
115                              warninp_t wi);
116 /* Even more checks, charge group radii vs. cut-off's only. */
117
118 void get_ir(const char *mdparin, const char *mdparout,
119             gmx::MDModules *mdModules, t_inputrec *ir, t_gromppopts *opts,
120             WriteMdpHeader writeMdpHeader, warninp_t wi);
121 /* Read the input file, and retrieve data for inputrec.
122  * More data are read, but the are only evaluated when the next
123  * function is called. Also prints the input file back to mdparout.
124  */
125
126 void do_index(const char* mdparin,
127               const char *ndx,
128               gmx_mtop_t *mtop,
129               bool        bVerbose,
130               t_inputrec *ir,
131               warninp_t   wi);
132 /* Read the index file and assign grp numbers to atoms.
133  */
134
135 /* Routines In readpull.c */
136
137 char **read_pullparams(std::vector<t_inpfile> *inp,
138                        pull_params_t          *pull,
139                        warninp_t               wi);
140 /* Reads the pull parameters, returns a list of the pull group names */
141
142 void make_pull_groups(pull_params_t *pull,
143                       char **pgnames,
144                       const t_blocka *grps, char **gnames);
145 /* Process the pull group parameters after reading the index groups */
146
147 void make_pull_coords(pull_params_t *pull);
148 /* Process the pull coordinates after reading the pull groups */
149
150 pull_t *set_pull_init(t_inputrec *ir, const gmx_mtop_t *mtop,
151                       rvec *x, matrix box, real lambda,
152                       warninp_t wi);
153 /* Prints the initial pull group distances in x.
154  * If requested, adds the current distance to the initial reference location.
155  * Returns the pull_t pull work struct. This should be passed to finish_pull()
156  * after all modules have registered their external potentials, if present.
157  */
158
159 char **read_rotparams(std::vector<t_inpfile> *inp, t_rot *rot, warninp_t wi);
160 /* Reads enforced rotation parameters, returns a list of the rot group names */
161
162 void make_rotation_groups(t_rot *rot, char **rotgnames,
163                           t_blocka *grps, char **gnames);
164 /* Process the rotation parameters after reading the index groups */
165
166 void set_reference_positions(t_rot *rot, rvec *x, matrix box,
167                              const char *fn, bool bSet, warninp_t wi);
168
169 #endif