2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
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"
59 static const char* RotStr = "Enforced rotation:";
62 static char s_vec[STRLEN];
65 static void string2dvec(char buf[], dvec nums)
67 if (sscanf(buf, "%lf%lf%lf", &nums[0], &nums[1], &nums[2]) != 3)
69 gmx_fatal(FARGS, "Expected three numbers at input line %s", buf);
74 extern std::vector<std::string> read_rotparams(std::vector<t_inpfile>* inp, t_rot* rot, warninp_t wi)
78 char warn_buf[STRLEN];
82 /* read rotation parameters */
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);
95 gmx_fatal(FARGS, "rot-ngroups should be >= 1");
98 snew(rot->grp, rot->ngrp);
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++)
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 = get_eenum(inp, buf, erotg_names);
116 printStringNoNewline(inp, "Use mass-weighting of the rotation group positions");
117 sprintf(buf, "rot-massw%d", g);
118 rotg->bMassW = get_eenum(inp, buf, yesno_names);
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 */
127 dsvmul(1.0 / dnorm(vec), vec, vec);
131 sprintf(warn_buf, "rot-vec%d = 0", g);
132 warning_error(wi, warn_buf);
135 "%s Group %d (%s) normalized rot. vector: %f %f %f\n",
138 erotg_names[rotg->eType],
142 for (m = 0; m < DIM; m++)
144 rotg->inputVec[m] = vec[m];
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");
151 if ((rotg->eType == erotgISO) || (rotg->eType == erotgPM) || (rotg->eType == erotgRM)
152 || (rotg->eType == erotgRM2))
154 string2dvec(s_vec, vec);
156 for (m = 0; m < DIM; m++)
158 rotg->pivot[m] = vec[m];
161 printStringNoNewline(inp, "Rotation rate (degree/ps) and force constant (kJ/(mol*nm^2))");
162 sprintf(buf, "rot-rate%d", g);
163 rotg->rate = get_ereal(inp, buf, 0.0, wi);
165 sprintf(buf, "rot-k%d", g);
166 rotg->k = get_ereal(inp, buf, 0.0, wi);
169 sprintf(warn_buf, "rot-k%d <= 0", g);
170 warning_note(wi, warn_buf);
173 printStringNoNewline(inp, "Slab distance for flexible axis rotation (nm)");
174 sprintf(buf, "rot-slab-dist%d", g);
175 rotg->slab_dist = get_ereal(inp, buf, 1.5, wi);
176 if (rotg->slab_dist <= 0.0)
178 sprintf(warn_buf, "rot-slab-dist%d <= 0", g);
179 warning_error(wi, warn_buf);
182 printStringNoNewline(inp,
183 "Minimum value of Gaussian function for the force to be evaluated "
184 "(for flex* potentials)");
185 sprintf(buf, "rot-min-gauss%d", g);
186 rotg->min_gaussian = get_ereal(inp, buf, 1e-3, wi);
187 if (rotg->min_gaussian <= 0.0)
189 sprintf(warn_buf, "rot-min-gauss%d <= 0", g);
190 warning_error(wi, warn_buf);
193 printStringNoNewline(
194 inp, "Value of additive constant epsilon' (nm^2) for rm2* and flex2* potentials");
195 sprintf(buf, "rot-eps%d", g);
196 rotg->eps = get_ereal(inp, buf, 1e-4, wi);
197 if ((rotg->eps <= 0.0) && (rotg->eType == erotgRM2 || rotg->eType == erotgFLEX2))
199 sprintf(warn_buf, "rot-eps%d <= 0", g);
200 warning_error(wi, warn_buf);
203 printStringNoNewline(
205 "Fitting method to determine angle of rotation group (rmsd, norm, or potential)");
206 sprintf(buf, "rot-fit-method%d", g);
207 rotg->eFittype = get_eenum(inp, buf, erotg_fitnames);
208 printStringNoNewline(inp,
209 "For fit type 'potential', nr. of angles around the reference for "
210 "which the pot. is evaluated");
211 sprintf(buf, "rot-potfit-nsteps%d", g);
212 rotg->PotAngle_nstep = get_eint(inp, buf, 21, wi);
213 if ((rotg->eFittype == erotgFitPOT) && (rotg->PotAngle_nstep < 1))
215 sprintf(warn_buf, "rot-potfit-nsteps%d < 1", g);
216 warning_error(wi, warn_buf);
218 printStringNoNewline(
219 inp, "For fit type 'potential', distance in degrees between two consecutive angles");
220 sprintf(buf, "rot-potfit-step%d", g);
221 rotg->PotAngle_step = get_ereal(inp, buf, 0.25, wi);
228 /* Check whether the box is unchanged */
229 static void check_box_unchanged(matrix f_box, matrix box, const char fn[], warninp_t wi)
233 char warn_buf[STRLEN];
236 for (i = 0; i < DIM; i++)
238 for (ii = 0; ii < DIM; ii++)
240 if (f_box[i][ii] != box[i][ii])
248 sprintf(warn_buf, "%s Box size in reference file %s differs from actual box size!", RotStr, fn);
249 warning(wi, warn_buf);
250 pr_rvecs(stderr, 0, "Your box is:", box, 3);
251 pr_rvecs(stderr, 0, "Box in file:", f_box, 3);
256 /* Extract the reference positions for the rotation group(s) */
257 extern void set_reference_positions(t_rot* rot, rvec* x, matrix box, const char* fn, bool bSet, warninp_t wi)
261 gmx_trr_header_t header; /* Header information of reference file */
262 rvec f_box[3]; /* Box from reference file */
264 for (g = 0; g < rot->ngrp; g++)
267 fprintf(stderr, "%s group %d has %d reference positions.\n", RotStr, g, rotg->nat);
268 snew(rotg->x_ref, rotg->nat);
270 /* Construct the name for the file containing the reference positions for this group: */
271 std::string reffileString =
272 gmx::Path::concatenateBeforeExtension(fn, gmx::formatString(".%d", g));
273 const char* reffile = reffileString.c_str();
275 /* If the base filename for the reference position files was explicitly set by
276 * the user, we issue a fatal error if the group file can not be found */
277 if (bSet && !gmx_fexist(reffile))
280 "%s The file containing the reference positions was not found.\n"
281 "Expected the file '%s' for group %d.\n",
287 if (gmx_fexist(reffile))
289 fprintf(stderr, " Reading them from %s.\n", reffile);
290 gmx_trr_read_single_header(reffile, &header);
291 if (rotg->nat != header.natoms)
294 "Number of atoms in file %s (%d) does not match the number of atoms in "
295 "rotation group (%d)!\n",
300 gmx_trr_read_single_frame(
301 reffile, &header.step, &header.t, &header.lambda, f_box, &header.natoms, rotg->x_ref, nullptr, nullptr);
303 /* Check whether the box is unchanged and output a warning if not: */
304 check_box_unchanged(f_box, box, reffile, wi);
308 fprintf(stderr, " Saving them to %s.\n", reffile);
309 for (i = 0; i < rotg->nat; i++)
312 copy_rvec(x[ii], rotg->x_ref[i]);
314 gmx_trr_write_single_frame(reffile, g, 0.0, 0.0, box, rotg->nat, rotg->x_ref, nullptr, nullptr);
320 extern void make_rotation_groups(t_rot* rot,
321 gmx::ArrayRef<const std::string> rotateGroupNames,
329 for (g = 0; g < rot->ngrp; g++)
332 ig = search_string(rotateGroupNames[g].c_str(), grps->nr, gnames);
333 rotg->nat = grps->index[ig + 1] - grps->index[ig];
337 fprintf(stderr, "Rotation group %d '%s' has %d atoms\n", g, rotateGroupNames[g].c_str(), rotg->nat);
338 snew(rotg->ind, rotg->nat);
339 for (i = 0; i < rotg->nat; i++)
341 rotg->ind[i] = grps->a[grps->index[ig] + i];
346 gmx_fatal(FARGS, "Rotation group %d '%s' is empty", g, rotateGroupNames[g].c_str());