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