Refactor md_enums
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / readrot.cpp
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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include <string>
41
42 #include "gromacs/fileio/readinp.h"
43 #include "gromacs/fileio/trrio.h"
44 #include "gromacs/fileio/warninp.h"
45 #include "gromacs/gmxpreprocess/readir.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/math/vecdump.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/topology/block.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/path.h"
56 #include "gromacs/utility/smalloc.h"
57 #include "gromacs/utility/stringutil.h"
58
59 static const char* RotStr = "Enforced rotation:";
60
61
62 static char s_vec[STRLEN];
63
64
65 static void string2dvec(char buf[], dvec nums)
66 {
67     if (sscanf(buf, "%lf%lf%lf", &nums[0], &nums[1], &nums[2]) != 3)
68     {
69         gmx_fatal(FARGS, "Expected three numbers at input line %s", buf);
70     }
71 }
72
73
74 extern std::vector<std::string> read_rotparams(std::vector<t_inpfile>* inp, t_rot* rot, warninp_t wi)
75 {
76     int       g, m;
77     char      buf[STRLEN];
78     char      warn_buf[STRLEN];
79     dvec      vec;
80     t_rotgrp* rotg;
81
82     /* read rotation parameters */
83     printStringNoNewline(
84             inp,
85             "Output frequency for angle, torque and rotation potential energy for the whole group");
86     rot->nstrout = get_eint(inp, "rot-nstrout", 100, wi);
87     printStringNoNewline(inp,
88                          "Output frequency for per-slab data (angles, torques and slab centers)");
89     rot->nstsout = get_eint(inp, "rot-nstsout", 1000, wi);
90     printStringNoNewline(inp, "Number of rotation groups");
91     rot->ngrp = get_eint(inp, "rot-ngroups", 1, wi);
92
93     if (rot->ngrp < 1)
94     {
95         gmx_fatal(FARGS, "rot-ngroups should be >= 1");
96     }
97
98     snew(rot->grp, rot->ngrp);
99
100     /* Read the rotation groups */
101     std::vector<std::string> rotateGroups(rot->ngrp);
102     char                     readBuffer[STRLEN];
103     for (g = 0; g < rot->ngrp; g++)
104     {
105         rotg = &rot->grp[g];
106         printStringNoNewline(inp, "Rotation group name");
107         sprintf(buf, "rot-group%d", g);
108         setStringEntry(inp, buf, readBuffer, "");
109         rotateGroups[g] = readBuffer;
110         printStringNoNewline(inp,
111                              "Rotation potential. Can be iso, iso-pf, pm, pm-pf, rm, rm-pf, rm2, "
112                              "rm2-pf, flex, flex-t, flex2, flex2-t");
113         sprintf(buf, "rot-type%d", g);
114         rotg->eType = getEnum<EnforcedRotationGroupType>(inp, buf, wi);
115
116         printStringNoNewline(inp, "Use mass-weighting of the rotation group positions");
117         sprintf(buf, "rot-massw%d", g);
118         rotg->bMassW = getEnum<Boolean>(inp, buf, wi) != Boolean::No;
119
120         printStringNoNewline(inp, "Rotation vector, will get normalized");
121         sprintf(buf, "rot-vec%d", g);
122         setStringEntry(inp, buf, s_vec, "1.0 0.0 0.0");
123         string2dvec(s_vec, vec);
124         /* Normalize the rotation vector */
125         if (dnorm(vec) != 0)
126         {
127             dsvmul(1.0 / dnorm(vec), vec, vec);
128         }
129         else
130         {
131             sprintf(warn_buf, "rot-vec%d = 0", g);
132             warning_error(wi, warn_buf);
133         }
134         fprintf(stderr,
135                 "%s Group %d (%s) normalized rot. vector: %f %f %f\n",
136                 RotStr,
137                 g,
138                 enumValueToString(rotg->eType),
139                 vec[0],
140                 vec[1],
141                 vec[2]);
142         for (m = 0; m < DIM; m++)
143         {
144             rotg->inputVec[m] = vec[m];
145         }
146
147         printStringNoNewline(inp, "Pivot point for the potentials iso, pm, rm, and rm2 (nm)");
148         sprintf(buf, "rot-pivot%d", g);
149         setStringEntry(inp, buf, s_vec, "0.0 0.0 0.0");
150         clear_dvec(vec);
151         if ((rotg->eType == EnforcedRotationGroupType::Iso)
152             || (rotg->eType == EnforcedRotationGroupType::Pm)
153             || (rotg->eType == EnforcedRotationGroupType::Rm)
154             || (rotg->eType == EnforcedRotationGroupType::Rm2))
155         {
156             string2dvec(s_vec, vec);
157         }
158         for (m = 0; m < DIM; m++)
159         {
160             rotg->pivot[m] = vec[m];
161         }
162
163         printStringNoNewline(inp, "Rotation rate (degree/ps) and force constant (kJ/(mol*nm^2))");
164         sprintf(buf, "rot-rate%d", g);
165         rotg->rate = get_ereal(inp, buf, 0.0, wi);
166
167         sprintf(buf, "rot-k%d", g);
168         rotg->k = get_ereal(inp, buf, 0.0, wi);
169         if (rotg->k <= 0.0)
170         {
171             sprintf(warn_buf, "rot-k%d <= 0", g);
172             warning_note(wi, warn_buf);
173         }
174
175         printStringNoNewline(inp, "Slab distance for flexible axis rotation (nm)");
176         sprintf(buf, "rot-slab-dist%d", g);
177         rotg->slab_dist = get_ereal(inp, buf, 1.5, wi);
178         if (rotg->slab_dist <= 0.0)
179         {
180             sprintf(warn_buf, "rot-slab-dist%d <= 0", g);
181             warning_error(wi, warn_buf);
182         }
183
184         printStringNoNewline(inp,
185                              "Minimum value of Gaussian function for the force to be evaluated "
186                              "(for flex* potentials)");
187         sprintf(buf, "rot-min-gauss%d", g);
188         rotg->min_gaussian = get_ereal(inp, buf, 1e-3, wi);
189         if (rotg->min_gaussian <= 0.0)
190         {
191             sprintf(warn_buf, "rot-min-gauss%d <= 0", g);
192             warning_error(wi, warn_buf);
193         }
194
195         printStringNoNewline(
196                 inp, "Value of additive constant epsilon' (nm^2) for rm2* and flex2* potentials");
197         sprintf(buf, "rot-eps%d", g);
198         rotg->eps = get_ereal(inp, buf, 1e-4, wi);
199         if ((rotg->eps <= 0.0)
200             && (rotg->eType == EnforcedRotationGroupType::Rm2
201                 || rotg->eType == EnforcedRotationGroupType::Flex2))
202         {
203             sprintf(warn_buf, "rot-eps%d <= 0", g);
204             warning_error(wi, warn_buf);
205         }
206
207         printStringNoNewline(
208                 inp,
209                 "Fitting method to determine angle of rotation group (rmsd, norm, or potential)");
210         sprintf(buf, "rot-fit-method%d", g);
211         rotg->eFittype = getEnum<RotationGroupFitting>(inp, buf, wi);
212         printStringNoNewline(inp,
213                              "For fit type 'potential', nr. of angles around the reference for "
214                              "which the pot. is evaluated");
215         sprintf(buf, "rot-potfit-nsteps%d", g);
216         rotg->PotAngle_nstep = get_eint(inp, buf, 21, wi);
217         if ((rotg->eFittype == RotationGroupFitting::Pot) && (rotg->PotAngle_nstep < 1))
218         {
219             sprintf(warn_buf, "rot-potfit-nsteps%d < 1", g);
220             warning_error(wi, warn_buf);
221         }
222         printStringNoNewline(
223                 inp, "For fit type 'potential', distance in degrees between two consecutive angles");
224         sprintf(buf, "rot-potfit-step%d", g);
225         rotg->PotAngle_step = get_ereal(inp, buf, 0.25, wi);
226     }
227
228     return rotateGroups;
229 }
230
231
232 /* Check whether the box is unchanged */
233 static void check_box_unchanged(matrix f_box, matrix box, const char fn[], warninp_t wi)
234 {
235     int  i, ii;
236     bool bSame = TRUE;
237     char warn_buf[STRLEN];
238
239
240     for (i = 0; i < DIM; i++)
241     {
242         for (ii = 0; ii < DIM; ii++)
243         {
244             if (f_box[i][ii] != box[i][ii])
245             {
246                 bSame = FALSE;
247             }
248         }
249     }
250     if (!bSame)
251     {
252         sprintf(warn_buf, "%s Box size in reference file %s differs from actual box size!", RotStr, fn);
253         warning(wi, warn_buf);
254         pr_rvecs(stderr, 0, "Your box is:", box, 3);
255         pr_rvecs(stderr, 0, "Box in file:", f_box, 3);
256     }
257 }
258
259
260 /* Extract the reference positions for the rotation group(s) */
261 extern void set_reference_positions(t_rot* rot, rvec* x, matrix box, const char* fn, bool bSet, warninp_t wi)
262 {
263     int              g, i, ii;
264     t_rotgrp*        rotg;
265     gmx_trr_header_t header;   /* Header information of reference file */
266     rvec             f_box[3]; /* Box from reference file */
267
268     for (g = 0; g < rot->ngrp; g++)
269     {
270         rotg = &rot->grp[g];
271         fprintf(stderr, "%s group %d has %d reference positions.\n", RotStr, g, rotg->nat);
272         snew(rotg->x_ref, rotg->nat);
273
274         /* Construct the name for the file containing the reference positions for this group: */
275         std::string reffileString =
276                 gmx::Path::concatenateBeforeExtension(fn, gmx::formatString(".%d", g));
277         const char* reffile = reffileString.c_str();
278
279         /* If the base filename for the reference position files was explicitly set by
280          * the user, we issue a fatal error if the group file can not be found */
281         if (bSet && !gmx_fexist(reffile))
282         {
283             gmx_fatal(FARGS,
284                       "%s The file containing the reference positions was not found.\n"
285                       "Expected the file '%s' for group %d.\n",
286                       RotStr,
287                       reffile,
288                       g);
289         }
290
291         if (gmx_fexist(reffile))
292         {
293             fprintf(stderr, "  Reading them from %s.\n", reffile);
294             gmx_trr_read_single_header(reffile, &header);
295             if (rotg->nat != header.natoms)
296             {
297                 gmx_fatal(FARGS,
298                           "Number of atoms in file %s (%d) does not match the number of atoms in "
299                           "rotation group (%d)!\n",
300                           reffile,
301                           header.natoms,
302                           rotg->nat);
303             }
304             gmx_trr_read_single_frame(
305                     reffile, &header.step, &header.t, &header.lambda, f_box, &header.natoms, rotg->x_ref, nullptr, nullptr);
306
307             /* Check whether the box is unchanged and output a warning if not: */
308             check_box_unchanged(f_box, box, reffile, wi);
309         }
310         else
311         {
312             fprintf(stderr, " Saving them to %s.\n", reffile);
313             for (i = 0; i < rotg->nat; i++)
314             {
315                 ii = rotg->ind[i];
316                 copy_rvec(x[ii], rotg->x_ref[i]);
317             }
318             gmx_trr_write_single_frame(reffile, g, 0.0, 0.0, box, rotg->nat, rotg->x_ref, nullptr, nullptr);
319         }
320     }
321 }
322
323
324 extern void make_rotation_groups(t_rot*                           rot,
325                                  gmx::ArrayRef<const std::string> rotateGroupNames,
326                                  t_blocka*                        grps,
327                                  char**                           gnames)
328 {
329     int       g, ig = -1, i;
330     t_rotgrp* rotg;
331
332
333     for (g = 0; g < rot->ngrp; g++)
334     {
335         rotg      = &rot->grp[g];
336         ig        = search_string(rotateGroupNames[g].c_str(), grps->nr, gnames);
337         rotg->nat = grps->index[ig + 1] - grps->index[ig];
338
339         if (rotg->nat > 0)
340         {
341             fprintf(stderr, "Rotation group %d '%s' has %d atoms\n", g, rotateGroupNames[g].c_str(), rotg->nat);
342             snew(rotg->ind, rotg->nat);
343             for (i = 0; i < rotg->nat; i++)
344             {
345                 rotg->ind[i] = grps->a[grps->index[ig] + i];
346             }
347         }
348         else
349         {
350             gmx_fatal(FARGS, "Rotation group %d '%s' is empty", g, rotateGroupNames[g].c_str());
351         }
352     }
353 }