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, 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.
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.
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.
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.
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.
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.
39 #include "gromacs/fileio/trnio.h"
40 #include "gromacs/gmxpreprocess/readir.h"
41 #include "gromacs/legacyheaders/names.h"
42 #include "gromacs/legacyheaders/txtdump.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/futil.h"
46 #include "gromacs/utility/smalloc.h"
48 static char *RotStr = {"Enforced rotation:"};
51 static char s_vec[STRLEN];
54 static void string2dvec(char buf[], dvec nums)
56 if (sscanf(buf, "%lf%lf%lf", &nums[0], &nums[1], &nums[2]) != 3)
58 gmx_fatal(FARGS, "Expected three numbers at input line %s", buf);
63 extern char **read_rotparams(int *ninp_p, t_inpfile **inp_p, t_rot *rot,
71 char warn_buf[STRLEN];
78 /* read rotation parameters */
79 CTYPE("Output frequency for angle, torque and rotation potential energy for the whole group");
80 ITYPE("rot_nstrout", rot->nstrout, 100);
81 CTYPE("Output frequency for per-slab data (angles, torques and slab centers)");
82 ITYPE("rot_nstsout", rot->nstsout, 1000);
83 CTYPE("Number of rotation groups");
84 ITYPE("rot_ngroups", rot->ngrp, 1);
88 gmx_fatal(FARGS, "rot_ngroups should be >= 1");
91 snew(rot->grp, rot->ngrp);
93 /* Read the rotation groups */
94 snew(grpbuf, rot->ngrp);
95 for (g = 0; g < rot->ngrp; g++)
98 snew(grpbuf[g], STRLEN);
99 CTYPE("Rotation group name");
100 sprintf(buf, "rot_group%d", g);
101 STYPE(buf, grpbuf[g], "");
103 CTYPE("Rotation potential. Can be iso, iso-pf, pm, pm-pf, rm, rm-pf, rm2, rm2-pf, flex, flex-t, flex2, flex2-t");
104 sprintf(buf, "rot_type%d", g);
105 ETYPE(buf, rotg->eType, erotg_names);
107 CTYPE("Use mass-weighting of the rotation group positions");
108 sprintf(buf, "rot_massw%d", g);
109 ETYPE(buf, rotg->bMassW, yesno_names);
111 CTYPE("Rotation vector, will get normalized");
112 sprintf(buf, "rot_vec%d", g);
113 STYPE(buf, s_vec, "1.0 0.0 0.0");
114 string2dvec(s_vec, vec);
115 /* Normalize the rotation vector */
118 dsvmul(1.0/dnorm(vec), vec, vec);
122 sprintf(warn_buf, "rot_vec%d = 0", g);
123 warning_error(wi, warn_buf);
125 fprintf(stderr, "%s Group %d (%s) normalized rot. vector: %f %f %f\n",
126 RotStr, g, erotg_names[rotg->eType], vec[0], vec[1], vec[2]);
127 for (m = 0; m < DIM; m++)
129 rotg->vec[m] = vec[m];
132 CTYPE("Pivot point for the potentials iso, pm, rm, and rm2 (nm)");
133 sprintf(buf, "rot_pivot%d", g);
134 STYPE(buf, s_vec, "0.0 0.0 0.0");
136 if ( (rotg->eType == erotgISO) || (rotg->eType == erotgPM) || (rotg->eType == erotgRM) || (rotg->eType == erotgRM2) )
138 string2dvec(s_vec, vec);
140 for (m = 0; m < DIM; m++)
142 rotg->pivot[m] = vec[m];
145 CTYPE("Rotation rate (degree/ps) and force constant (kJ/(mol*nm^2))");
146 sprintf(buf, "rot_rate%d", g);
147 RTYPE(buf, rotg->rate, 0.0);
149 sprintf(buf, "rot_k%d", g);
150 RTYPE(buf, rotg->k, 0.0);
153 sprintf(warn_buf, "rot_k%d <= 0", g);
154 warning_note(wi, warn_buf);
157 CTYPE("Slab distance for flexible axis rotation (nm)");
158 sprintf(buf, "rot_slab_dist%d", g);
159 RTYPE(buf, rotg->slab_dist, 1.5);
160 if (rotg->slab_dist <= 0.0)
162 sprintf(warn_buf, "rot_slab_dist%d <= 0", g);
163 warning_error(wi, warn_buf);
166 CTYPE("Minimum value of Gaussian function for the force to be evaluated (for flex* potentials)");
167 sprintf(buf, "rot_min_gauss%d", g);
168 RTYPE(buf, rotg->min_gaussian, 1e-3);
169 if (rotg->min_gaussian <= 0.0)
171 sprintf(warn_buf, "rot_min_gauss%d <= 0", g);
172 warning_error(wi, warn_buf);
175 CTYPE("Value of additive constant epsilon' (nm^2) for rm2* and flex2* potentials");
176 sprintf(buf, "rot_eps%d", g);
177 RTYPE(buf, rotg->eps, 1e-4);
178 if ( (rotg->eps <= 0.0) && (rotg->eType == erotgRM2 || rotg->eType == erotgFLEX2) )
180 sprintf(warn_buf, "rot_eps%d <= 0", g);
181 warning_error(wi, warn_buf);
184 CTYPE("Fitting method to determine angle of rotation group (rmsd, norm, or potential)");
185 sprintf(buf, "rot_fit_method%d", g);
186 ETYPE(buf, rotg->eFittype, erotg_fitnames);
187 CTYPE("For fit type 'potential', nr. of angles around the reference for which the pot. is evaluated");
188 sprintf(buf, "rot_potfit_nsteps%d", g);
189 ITYPE(buf, rotg->PotAngle_nstep, 21);
190 if ( (rotg->eFittype == erotgFitPOT) && (rotg->PotAngle_nstep < 1) )
192 sprintf(warn_buf, "rot_potfit_nsteps%d < 1", g);
193 warning_error(wi, warn_buf);
195 CTYPE("For fit type 'potential', distance in degrees between two consecutive angles");
196 sprintf(buf, "rot_potfit_step%d", g);
197 RTYPE(buf, rotg->PotAngle_step, 0.25);
207 /* Check whether the box is unchanged */
208 static void check_box_unchanged(matrix f_box, matrix box, char fn[], warninp_t wi)
211 gmx_bool bSame = TRUE;
212 char warn_buf[STRLEN];
215 for (i = 0; i < DIM; i++)
217 for (ii = 0; ii < DIM; ii++)
219 if (f_box[i][ii] != box[i][ii])
227 sprintf(warn_buf, "%s Box size in reference file %s differs from actual box size!",
229 warning(wi, warn_buf);
230 pr_rvecs(stderr, 0, "Your box is:", box, 3);
231 pr_rvecs(stderr, 0, "Box in file:", f_box, 3);
236 /* Extract the reference positions for the rotation group(s) */
237 extern void set_reference_positions(
238 t_rot *rot, rvec *x, matrix box,
239 const char *fn, gmx_bool bSet, warninp_t wi)
243 t_trnheader header; /* Header information of reference file */
244 char base[STRLEN], extension[STRLEN], reffile[STRLEN];
246 rvec f_box[3]; /* Box from reference file */
249 /* Base name and extension of the reference file: */
250 strncpy(base, fn, STRLEN - 1);
251 base[STRLEN-1] = '\0';
252 extpos = strrchr(base, '.');
253 strcpy(extension, extpos+1);
257 for (g = 0; g < rot->ngrp; g++)
260 fprintf(stderr, "%s group %d has %d reference positions.\n", RotStr, g, rotg->nat);
261 snew(rotg->x_ref, rotg->nat);
263 /* Construct the name for the file containing the reference positions for this group: */
264 sprintf(reffile, "%s.%d.%s", base, g, extension);
266 /* If the base filename for the reference position files was explicitly set by
267 * the user, we issue a fatal error if the group file can not be found */
268 if (bSet && !gmx_fexist(reffile))
270 gmx_fatal(FARGS, "%s The file containing the reference positions was not found.\n"
271 "Expected the file '%s' for group %d.\n",
275 if (gmx_fexist(reffile))
277 fprintf(stderr, " Reading them from %s.\n", reffile);
278 read_trnheader(reffile, &header);
279 if (rotg->nat != header.natoms)
281 gmx_fatal(FARGS, "Number of atoms in file %s (%d) does not match the number of atoms in rotation group (%d)!\n",
282 reffile, header.natoms, rotg->nat);
284 read_trn(reffile, &header.step, &header.t, &header.lambda, f_box, &header.natoms, rotg->x_ref, NULL, NULL);
286 /* Check whether the box is unchanged and output a warning if not: */
287 check_box_unchanged(f_box, box, reffile, wi);
291 fprintf(stderr, " Saving them to %s.\n", reffile);
292 for (i = 0; i < rotg->nat; i++)
295 copy_rvec(x[ii], rotg->x_ref[i]);
297 write_trn(reffile, g, 0.0, 0.0, box, rotg->nat, rotg->x_ref, NULL, NULL);
303 extern void make_rotation_groups(t_rot *rot, char **rotgnames, t_blocka *grps, char **gnames)
309 for (g = 0; g < rot->ngrp; g++)
312 ig = search_string(rotgnames[g], grps->nr, gnames);
313 rotg->nat = grps->index[ig+1] - grps->index[ig];
317 fprintf(stderr, "Rotation group %d '%s' has %d atoms\n", g, rotgnames[g], rotg->nat);
318 snew(rotg->ind, rotg->nat);
319 for (i = 0; i < rotg->nat; i++)
321 rotg->ind[i] = grps->a[grps->index[ig]+i];
326 gmx_fatal(FARGS, "Rotation group %d '%s' is empty", g, rotgnames[g]);