70cfc5e4c9c1e8749b2cfb427dbdd5e7de5dbc5b
[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,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 #include "gmxpre.h"
38
39 #include "gromacs/fileio/readinp.h"
40 #include "gromacs/fileio/trrio.h"
41 #include "gromacs/fileio/warninp.h"
42 #include "gromacs/gmxpreprocess/readir.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/math/vecdump.h"
45 #include "gromacs/mdtypes/inputrec.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/topology/block.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/futil.h"
51 #include "gromacs/utility/smalloc.h"
52
53 static const char *RotStr = "Enforced rotation:";
54
55
56 static char s_vec[STRLEN];
57
58
59 static void string2dvec(char buf[], dvec nums)
60 {
61     if (sscanf(buf, "%lf%lf%lf", &nums[0], &nums[1], &nums[2]) != 3)
62     {
63         gmx_fatal(FARGS, "Expected three numbers at input line %s", buf);
64     }
65 }
66
67
68 extern char **read_rotparams(std::vector<t_inpfile> *inp, t_rot *rot,
69                              warninp_t wi)
70 {
71     int                    g, m;
72     char                 **grpbuf;
73     char                   buf[STRLEN];
74     char                   warn_buf[STRLEN];
75     dvec                   vec;
76     t_rotgrp              *rotg;
77
78     /* read rotation parameters */
79     printStringNoNewline(inp, "Output frequency for angle, torque and rotation potential energy for the whole group");
80     rot->nstrout = get_eint(inp, "rot-nstrout", 100, wi);
81     printStringNoNewline(inp, "Output frequency for per-slab data (angles, torques and slab centers)");
82     rot->nstsout = get_eint(inp, "rot-nstsout", 1000, wi);
83     printStringNoNewline(inp, "Number of rotation groups");
84     rot->ngrp = get_eint(inp, "rot-ngroups", 1, wi);
85
86     if (rot->ngrp < 1)
87     {
88         gmx_fatal(FARGS, "rot-ngroups should be >= 1");
89     }
90
91     snew(rot->grp, rot->ngrp);
92
93     /* Read the rotation groups */
94     snew(grpbuf, rot->ngrp);
95     for (g = 0; g < rot->ngrp; g++)
96     {
97         rotg = &rot->grp[g];
98         snew(grpbuf[g], STRLEN);
99         printStringNoNewline(inp, "Rotation group name");
100         sprintf(buf, "rot-group%d", g);
101         setStringEntry(inp, buf, grpbuf[g], "");
102
103         printStringNoNewline(inp, "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         rotg->eType = get_eenum(inp, buf, erotg_names);
106
107         printStringNoNewline(inp, "Use mass-weighting of the rotation group positions");
108         sprintf(buf, "rot-massw%d", g);
109         rotg->bMassW = get_eenum(inp, buf, yesno_names);
110
111         printStringNoNewline(inp, "Rotation vector, will get normalized");
112         sprintf(buf, "rot-vec%d", g);
113         setStringEntry(inp, buf, s_vec, "1.0 0.0 0.0");
114         string2dvec(s_vec, vec);
115         /* Normalize the rotation vector */
116         if (dnorm(vec) != 0)
117         {
118             dsvmul(1.0/dnorm(vec), vec, vec);
119         }
120         else
121         {
122             sprintf(warn_buf, "rot-vec%d = 0", g);
123             warning_error(wi, warn_buf);
124         }
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++)
128         {
129             rotg->vec[m] = vec[m];
130         }
131
132         printStringNoNewline(inp, "Pivot point for the potentials iso, pm, rm, and rm2 (nm)");
133         sprintf(buf, "rot-pivot%d", g);
134         setStringEntry(inp, buf, s_vec, "0.0 0.0 0.0");
135         clear_dvec(vec);
136         if ( (rotg->eType == erotgISO) || (rotg->eType == erotgPM) || (rotg->eType == erotgRM) || (rotg->eType == erotgRM2) )
137         {
138             string2dvec(s_vec, vec);
139         }
140         for (m = 0; m < DIM; m++)
141         {
142             rotg->pivot[m] = vec[m];
143         }
144
145         printStringNoNewline(inp, "Rotation rate (degree/ps) and force constant (kJ/(mol*nm^2))");
146         sprintf(buf, "rot-rate%d", g);
147         rotg->rate = get_ereal(inp, buf, 0.0, wi);
148
149         sprintf(buf, "rot-k%d", g);
150         rotg->k = get_ereal(inp, buf, 0.0, wi);
151         if (rotg->k <= 0.0)
152         {
153             sprintf(warn_buf, "rot-k%d <= 0", g);
154             warning_note(wi, warn_buf);
155         }
156
157         printStringNoNewline(inp, "Slab distance for flexible axis rotation (nm)");
158         sprintf(buf, "rot-slab-dist%d", g);
159         rotg->slab_dist = get_ereal(inp, buf, 1.5, wi);
160         if (rotg->slab_dist <= 0.0)
161         {
162             sprintf(warn_buf, "rot-slab-dist%d <= 0", g);
163             warning_error(wi, warn_buf);
164         }
165
166         printStringNoNewline(inp, "Minimum value of Gaussian function for the force to be evaluated (for flex* potentials)");
167         sprintf(buf, "rot-min-gauss%d", g);
168         rotg->min_gaussian = get_ereal(inp, buf, 1e-3, wi);
169         if (rotg->min_gaussian <= 0.0)
170         {
171             sprintf(warn_buf, "rot-min-gauss%d <= 0", g);
172             warning_error(wi, warn_buf);
173         }
174
175         printStringNoNewline(inp, "Value of additive constant epsilon' (nm^2) for rm2* and flex2* potentials");
176         sprintf(buf, "rot-eps%d", g);
177         rotg->eps = get_ereal(inp, buf, 1e-4, wi);
178         if ( (rotg->eps <= 0.0) && (rotg->eType == erotgRM2 || rotg->eType == erotgFLEX2) )
179         {
180             sprintf(warn_buf, "rot-eps%d <= 0", g);
181             warning_error(wi, warn_buf);
182         }
183
184         printStringNoNewline(inp, "Fitting method to determine angle of rotation group (rmsd, norm, or potential)");
185         sprintf(buf, "rot-fit-method%d", g);
186         rotg->eFittype = get_eenum(inp, buf, erotg_fitnames);
187         printStringNoNewline(inp, "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         rotg->PotAngle_nstep = get_eint(inp, buf, 21, wi);
190         if ( (rotg->eFittype == erotgFitPOT) && (rotg->PotAngle_nstep < 1) )
191         {
192             sprintf(warn_buf, "rot-potfit-nsteps%d < 1", g);
193             warning_error(wi, warn_buf);
194         }
195         printStringNoNewline(inp, "For fit type 'potential', distance in degrees between two consecutive angles");
196         sprintf(buf, "rot-potfit-step%d", g);
197         rotg->PotAngle_step = get_ereal(inp, buf, 0.25, wi);
198     }
199
200     return grpbuf;
201 }
202
203
204 /* Check whether the box is unchanged */
205 static void check_box_unchanged(matrix f_box, matrix box, char fn[], warninp_t wi)
206 {
207     int      i, ii;
208     gmx_bool bSame = TRUE;
209     char     warn_buf[STRLEN];
210
211
212     for (i = 0; i < DIM; i++)
213     {
214         for (ii = 0; ii < DIM; ii++)
215         {
216             if (f_box[i][ii] != box[i][ii])
217             {
218                 bSame = FALSE;
219             }
220         }
221     }
222     if (!bSame)
223     {
224         sprintf(warn_buf, "%s Box size in reference file %s differs from actual box size!",
225                 RotStr, fn);
226         warning(wi, warn_buf);
227         pr_rvecs(stderr, 0, "Your box is:", box, 3);
228         pr_rvecs(stderr, 0, "Box in file:", f_box, 3);
229     }
230 }
231
232
233 /* Extract the reference positions for the rotation group(s) */
234 extern void set_reference_positions(
235         t_rot *rot, rvec *x, matrix box,
236         const char *fn, gmx_bool bSet, warninp_t wi)
237 {
238     int              g, i, ii;
239     t_rotgrp        *rotg;
240     gmx_trr_header_t header;    /* Header information of reference file */
241     char             base[STRLEN], extension[STRLEN], reffile[STRLEN];
242     char            *extpos;
243     rvec             f_box[3]; /* Box from reference file */
244
245
246     /* Base name and extension of the reference file: */
247     strncpy(base, fn, STRLEN - 1);
248     base[STRLEN-1] = '\0';
249     extpos         = strrchr(base, '.');
250     strcpy(extension, extpos+1);
251     *extpos = '\0';
252
253
254     for (g = 0; g < rot->ngrp; g++)
255     {
256         rotg = &rot->grp[g];
257         fprintf(stderr, "%s group %d has %d reference positions.\n", RotStr, g, rotg->nat);
258         snew(rotg->x_ref, rotg->nat);
259
260         /* Construct the name for the file containing the reference positions for this group: */
261         sprintf(reffile, "%s.%d.%s", base, g, extension);
262
263         /* If the base filename for the reference position files was explicitly set by
264          * the user, we issue a fatal error if the group file can not be found */
265         if (bSet && !gmx_fexist(reffile))
266         {
267             gmx_fatal(FARGS, "%s The file containing the reference positions was not found.\n"
268                       "Expected the file '%s' for group %d.\n",
269                       RotStr, reffile, g);
270         }
271
272         if (gmx_fexist(reffile))
273         {
274             fprintf(stderr, "  Reading them from %s.\n", reffile);
275             gmx_trr_read_single_header(reffile, &header);
276             if (rotg->nat != header.natoms)
277             {
278                 gmx_fatal(FARGS, "Number of atoms in file %s (%d) does not match the number of atoms in rotation group (%d)!\n",
279                           reffile, header.natoms, rotg->nat);
280             }
281             gmx_trr_read_single_frame(reffile, &header.step, &header.t, &header.lambda, f_box, &header.natoms, rotg->x_ref, nullptr, nullptr);
282
283             /* Check whether the box is unchanged and output a warning if not: */
284             check_box_unchanged(f_box, box, reffile, wi);
285         }
286         else
287         {
288             fprintf(stderr, " Saving them to %s.\n", reffile);
289             for (i = 0; i < rotg->nat; i++)
290             {
291                 ii = rotg->ind[i];
292                 copy_rvec(x[ii], rotg->x_ref[i]);
293             }
294             gmx_trr_write_single_frame(reffile, g, 0.0, 0.0, box, rotg->nat, rotg->x_ref, nullptr, nullptr);
295         }
296     }
297 }
298
299
300 extern void make_rotation_groups(t_rot *rot, char **rotgnames, t_blocka *grps, char **gnames)
301 {
302     int       g, ig = -1, i;
303     t_rotgrp *rotg;
304
305
306     for (g = 0; g < rot->ngrp; g++)
307     {
308         rotg      = &rot->grp[g];
309         ig        = search_string(rotgnames[g], grps->nr, gnames);
310         rotg->nat = grps->index[ig+1] - grps->index[ig];
311
312         if (rotg->nat > 0)
313         {
314             fprintf(stderr, "Rotation group %d '%s' has %d atoms\n", g, rotgnames[g], rotg->nat);
315             snew(rotg->ind, rotg->nat);
316             for (i = 0; i < rotg->nat; i++)
317             {
318                 rotg->ind[i] = grps->a[grps->index[ig]+i];
319             }
320         }
321         else
322         {
323             gmx_fatal(FARGS, "Rotation group %d '%s' is empty", g, rotgnames[g]);
324         }
325     }
326 }