4a8e521762fe7d2ee9eab6abd36f14793be6a43d
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_make_edi.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /* The make_edi program was generously contributed by Oliver Lange, based
36  * on the code from g_anaeig. You can reach him as olange@gwdg.de. He
37  * probably also holds copyright to the following code.
38  */
39 #include "gmxpre.h"
40
41 #include <cctype>
42 #include <cmath>
43 #include <cstdlib>
44 #include <cstring>
45
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/pdbio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/gmxana/eigio.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61
62 typedef struct
63 {
64     real        deltaF0;
65     gmx_bool    bHarmonic;
66     gmx_bool    bConstForce;   /* Do constant force flooding instead of
67                                   evaluating a flooding potential             */
68     real        tau;
69     real        deltaF;
70     real        kT;
71     real        constEfl;
72     real        alpha2;
73 } t_edflood;
74
75
76 /* This type is for the average, reference, target, and origin structure   */
77 typedef struct edix
78 {
79     int          nr;            /* number of atoms this structure contains */
80     int         *anrs;          /* atom index numbers                      */
81     rvec        *x;             /* positions                               */
82     real        *sqrtm;         /* sqrt of the masses used for mass-
83                                  * weighting of analysis                   */
84 } t_edix;
85
86
87 typedef struct edipar
88 {
89     int         nini;           /* total Nr of atoms                    */
90     gmx_bool    fitmas;         /* true if trans fit with cm            */
91     gmx_bool    pcamas;         /* true if mass-weighted PCA            */
92     int         presteps;       /* number of steps to run without any
93                                  *    perturbations ... just monitoring */
94     int         outfrq;         /* freq (in steps) of writing to edo    */
95     int         maxedsteps;     /* max nr of steps per cycle            */
96     struct edix sref;           /* reference positions, to these fitting
97                                  * will be done                         */
98     struct edix sav;            /* average positions                    */
99     struct edix star;           /* target positions                     */
100     struct edix sori;           /* origin positions                     */
101     real        slope;          /* minimal slope in acceptance radexp   */
102     int         ned;            /* Nr of atoms in essdyn buffer         */
103     t_edflood   flood;          /* parameters especially for flooding   */
104 } t_edipar;
105
106
107
108 static void make_t_edx(struct edix *edx, int natoms, rvec *pos, int index[])
109 {
110     edx->nr   = natoms;
111     edx->anrs = index;
112     edx->x    = pos;
113 }
114
115 static void write_t_edx(FILE *fp, struct edix edx, const char *comment)
116 {
117     /*here we copy only the pointers into the t_edx struct
118        no data is copied and edx.box is ignored  */
119     int i;
120     fprintf(fp, "#%s \n %d \n", comment, edx.nr);
121     for (i = 0; i < edx.nr; i++)
122     {
123         fprintf(fp, "%d  %f  %f  %f\n", (edx.anrs)[i]+1, (edx.x)[i][XX], (edx.x)[i][YY], (edx.x)[i][ZZ]);
124     }
125 }
126
127 static int sscan_list(int *list[], const char *str, const char *listname)
128 {
129     /*this routine scans a string of the form 1,3-6,9 and returns the
130        selected numbers (in this case 1 3 4 5 6 9) in NULL-terminated array of integers.
131        memory for this list will be allocated  in this routine -- sscan_list expects *list to
132        be a NULL-Pointer
133
134        listname is a string used in the errormessage*/
135
136
137     int   i, istep;
138     char  c;
139     char *pos, *startpos, *step;
140     int   n = std::strlen(str);
141
142     /*enums to define the different lexical stati */
143     enum {
144         sBefore, sNumber, sMinus, sRange, sZero, sSmaller, sError, sSteppedRange
145     };
146
147     int   status     = sBefore; /*status of the deterministic automat to scan str   */
148     int   number     = 0;
149     int   end_number = 0;
150
151     char *start = nullptr; /*holds the string of the number behind a ','*/
152     char *end   = nullptr; /*holds the string of the number behind a '-' */
153
154     int   nvecs = 0;       /* counts the number of vectors in the list*/
155
156     step = nullptr;
157     snew(pos, n+4);
158     startpos = pos;
159     std::strcpy(pos, str);
160     pos[n]   = ',';
161     pos[n+1] = '1';
162     pos[n+2] = '\0';
163
164     *list = nullptr;
165
166     while ((c = *pos) != 0)
167     {
168         switch (status)
169         {
170             /* expect a number */
171             case sBefore:
172                 if (std::isdigit(c))
173                 {
174                     start  = pos;
175                     status = sNumber;
176                     break;
177                 }
178                 else
179                 {
180                     status = sError;
181                 }
182                 break;
183
184             /* have read a number, expect ',' or '-' */
185             case sNumber:
186                 if (c == ',')
187                 {
188                     /*store number*/
189                     srenew(*list, nvecs+1);
190                     (*list)[nvecs++] = number = std::strtol(start, nullptr, 10);
191                     status           = sBefore;
192                     if (number == 0)
193                     {
194                         status = sZero;
195                     }
196                     break;
197                 }
198                 else if (c == '-')
199                 {
200                     status = sMinus; break;
201                 }
202                 else if (std::isdigit(c))
203                 {
204                     break;
205                 }
206                 else
207                 {
208                     status = sError;
209                 }
210                 break;
211
212             /* have read a '-' -> expect a number */
213             case sMinus:
214                 if (std::isdigit(c))
215                 {
216                     end    = pos;
217                     status = sRange; break;
218                 }
219                 else
220                 {
221                     status = sError;
222                 }
223                 break;
224
225             case sSteppedRange:
226                 if (std::isdigit(c))
227                 {
228                     if (step)
229                     {
230                         status = sError; break;
231                     }
232                     else
233                     {
234                         step = pos;
235                     }
236                     status = sRange;
237                     break;
238                 }
239                 else
240                 {
241                     status = sError;
242                 }
243                 break;
244
245             /* have read the number after a minus, expect ',' or ':' */
246             case sRange:
247                 if (c == ',')
248                 {
249                     /*store numbers*/
250                     end_number = std::strtol(end, nullptr, 10);
251                     number     = std::strtol(start, nullptr, 10);
252                     status     = sBefore;
253                     if (number == 0)
254                     {
255                         status = sZero; break;
256                     }
257                     if (end_number <= number)
258                     {
259                         status = sSmaller; break;
260                     }
261                     srenew(*list, nvecs+end_number-number+1);
262                     if (step)
263                     {
264                         istep = strtol(step, nullptr, 10);
265                         step  = nullptr;
266                     }
267                     else
268                     {
269                         istep = 1;
270                     }
271                     for (i = number; i <= end_number; i += istep)
272                     {
273                         (*list)[nvecs++] = i;
274                     }
275                     break;
276                 }
277                 else if (c == ':')
278                 {
279                     status = sSteppedRange;
280                     break;
281                 }
282                 else if (std::isdigit(c))
283                 {
284                     break;
285                 }
286                 else
287                 {
288                     status = sError;
289                 }
290                 break;
291
292             /* format error occured */
293             case sError:
294                 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %d with char %c", listname, pos-startpos, *(pos-1));
295                 break;
296             /* logical error occured */
297             case sZero:
298                 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %d: eigenvector 0 is not valid", listname, pos-startpos);
299                 break;
300             case sSmaller:
301                 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %d: second index %d is not bigger than %d", listname, pos-startpos, end_number, number);
302                 break;
303         }
304         ++pos; /* read next character */
305     }          /*scanner has finished */
306
307     /* append zero to list of eigenvectors */
308     srenew(*list, nvecs+1);
309     (*list)[nvecs] = 0;
310     sfree(startpos);
311     return nvecs;
312 } /*sscan_list*/
313
314 static void write_eigvec(FILE* fp, int natoms, int eig_list[], rvec** eigvecs, int nvec, const char *grouptitle, real steps[])
315 {
316 /* eig_list is a zero-terminated list of indices into the eigvecs array.
317    eigvecs are coordinates of eigenvectors
318    grouptitle to write in the comment line
319    steps  -- array with stepsizes for evLINFIX, evLINACC and evRADACC
320  */
321
322     int  n = 0, i; rvec x;
323     while (eig_list[n++])
324     {
325         ;                 /*count selected eigenvecs*/
326
327     }
328     fprintf(fp, "# NUMBER OF EIGENVECTORS + %s\n %d\n", grouptitle, n-1);
329
330     /* write list of eigenvector indicess */
331     for (n = 0; eig_list[n]; n++)
332     {
333         if (steps)
334         {
335             fprintf(fp, "%8d   %g\n", eig_list[n], steps[n]);
336         }
337         else
338         {
339             fprintf(fp, "%8d   %g\n", eig_list[n], 1.0);
340         }
341     }
342     n = 0;
343
344     /* dump coordinates of the selected eigenvectors */
345     while (eig_list[n])
346     {
347         for (i = 0; i < natoms; i++)
348         {
349             if (eig_list[n] > nvec)
350             {
351                 gmx_fatal(FARGS, "Selected eigenvector %d is higher than maximum number %d of available eigenvectors", eig_list[n], nvec);
352             }
353             copy_rvec(eigvecs[eig_list[n]-1][i], x);
354             fprintf(fp, "%8.5f %8.5f %8.5f\n", x[XX], x[YY], x[ZZ]);
355         }
356         n++;
357     }
358 }
359
360
361 /*enum referring to the different lists of eigenvectors*/
362 enum {
363     evLINFIX, evLINACC, evFLOOD, evRADFIX, evRADACC, evRADCON, evMON,  evNr
364 };
365 #define oldMAGIC 666
366 #define MAGIC 670
367
368
369 static void write_the_whole_thing(FILE* fp, t_edipar *edpars, rvec** eigvecs,
370                                   int nvec, int *eig_listen[], real* evStepList[])
371 {
372 /* write edi-file */
373
374     /*Header*/
375     fprintf(fp, "#MAGIC\n %d \n#NINI\n %d\n#FITMAS\n %d\n#ANALYSIS_MAS\n %d\n",
376             MAGIC, edpars->nini, edpars->fitmas, edpars->pcamas);
377     fprintf(fp, "#OUTFRQ\n %d\n#MAXLEN\n %d\n#SLOPECRIT\n %f\n",
378             edpars->outfrq, edpars->maxedsteps, edpars->slope);
379     fprintf(fp, "#PRESTEPS\n %d\n#DELTA_F0\n %f\n#INIT_DELTA_F\n %f\n#TAU\n %f\n#EFL_NULL\n %f\n#ALPHA2\n %f\n#KT\n %f\n#HARMONIC\n %d\n#CONST_FORCE_FLOODING\n %d\n",
380             edpars->presteps, edpars->flood.deltaF0, edpars->flood.deltaF, edpars->flood.tau, edpars->flood.constEfl,
381             edpars->flood.alpha2, edpars->flood.kT, edpars->flood.bHarmonic, edpars->flood.bConstForce);
382
383     /* Average and reference positions */
384     write_t_edx(fp, edpars->sref, "NREF, XREF");
385     write_t_edx(fp, edpars->sav, "NAV, XAV");
386
387     /*Eigenvectors */
388
389     write_eigvec(fp, edpars->ned, eig_listen[evMON], eigvecs, nvec, "COMPONENTS GROUP 1", nullptr);
390     write_eigvec(fp, edpars->ned, eig_listen[evLINFIX], eigvecs, nvec, "COMPONENTS GROUP 2", evStepList[evLINFIX]);
391     write_eigvec(fp, edpars->ned, eig_listen[evLINACC], eigvecs, nvec, "COMPONENTS GROUP 3", evStepList[evLINACC]);
392     write_eigvec(fp, edpars->ned, eig_listen[evRADFIX], eigvecs, nvec, "COMPONENTS GROUP 4", evStepList[evRADFIX]);
393     write_eigvec(fp, edpars->ned, eig_listen[evRADACC], eigvecs, nvec, "COMPONENTS GROUP 5", nullptr);
394     write_eigvec(fp, edpars->ned, eig_listen[evRADCON], eigvecs, nvec, "COMPONENTS GROUP 6", nullptr);
395     write_eigvec(fp, edpars->ned, eig_listen[evFLOOD], eigvecs, nvec, "COMPONENTS GROUP 7", evStepList[evFLOOD]);
396
397
398     /*Target and Origin positions */
399     write_t_edx(fp, edpars->star, "NTARGET, XTARGET");
400     write_t_edx(fp, edpars->sori, "NORIGIN, XORIGIN");
401 }
402
403 static int read_conffile(const char *confin, rvec **x)
404 {
405     t_topology  top;
406     matrix      box;
407     printf("read coordnumber from file %s\n", confin);
408     read_tps_conf(confin, &top, nullptr, x, nullptr, box, FALSE);
409     printf("number of coordinates in file %d\n", top.atoms.nr);
410     return top.atoms.nr;
411 }
412
413
414 static void read_eigenvalues(int vecs[], const char *eigfile, real values[],
415                              gmx_bool bHesse, real kT, int natoms_average_struct)
416 {
417     int      neig, nrow, i;
418     double **eigval;
419
420     neig = read_xvg(eigfile, &eigval, &nrow);
421
422     fprintf(stderr, "Read %d eigenvalues\n", neig);
423     for (i = bHesse ? 6 : 0; i < neig; i++)
424     {
425         if (eigval[1][i] < -0.001 && bHesse)
426         {
427             fprintf(stderr,
428                     "WARNING: The Hessian Matrix has negative eigenvalue %f, we set it to zero (no flooding in this direction)\n\n", eigval[1][i]);
429         }
430
431         if (eigval[1][i] < 0)
432         {
433             eigval[1][i] = 0;
434         }
435     }
436     if (bHesse)
437     {
438         for (i = 0; vecs[i]; i++)
439         {
440             if (vecs[i] < 7)
441             {
442                 gmx_fatal(FARGS, "ERROR: You have chosen one of the first 6 eigenvectors of the HESSE Matrix. That does not make sense, since they correspond to the 6 rotational and translational degrees of freedom.\n\n");
443             }
444             values[i] = eigval[1][vecs[i]-1]/kT;
445         }
446     }
447     else
448     {
449         for (i = 0; vecs[i]; i++)
450         {
451             /* Make sure this eigenvalue does not correspond to one of the last 6 eigenvectors of the
452              * covariance matrix. These correspond to the rotational and translational degrees of
453              * freedom and will be zero within numerical accuracy.
454              *
455              * Note that the total number of eigenvectors produced by gmx covar depends on:
456              * 1) the total number of degrees of freedom of the system (3N, with N the number of atoms)
457              * 2) the number S of independent configurations fed into gmx covar.
458              * For long trajectories with lots of frames, usually S >= 3N + 1, so that one indeed gets
459              * 3N eigenvalues (of which the last 6 will have zero eigenvalues).
460              * For S < 3N + 1, however, the covariance matrix becomes rank deficient, and the number
461              * of possible eigenvalues is just S - 1. Since in make_edi we only know N but not S, we can
462              * only warn the user if he picked one of the last 6 of 3N.
463              */
464             if (vecs[i] > ( 3 * natoms_average_struct - 6 ))
465             {
466                 gmx_fatal(FARGS, "ERROR: You have chosen one of the last 6 eigenvectors of the COVARIANCE Matrix. That does not make sense, since they correspond to the 6 rotational and translational degrees of freedom.\n\n");
467             }
468             values[i] = 1/eigval[1][vecs[i]-1];
469         }
470     }
471     /* free memory */
472     for (i = 0; i < nrow; i++)
473     {
474         sfree(eigval[i]);
475     }
476     sfree(eigval);
477 }
478
479
480 static real *scan_vecparams(const char *str, const char * par, int nvecs)
481 {
482     char    f0[256], f1[256];         /*format strings adapted every pass of the loop*/
483     double  d;
484     int     i;
485     real   *vec_params;
486
487     snew(vec_params, nvecs);
488     if (str)
489     {
490         f0[0] = '\0';
491         for (i = 0; (i < nvecs); i++)
492         {
493             std::strcpy(f1, f0);    /*f0 is the format string for the "to-be-ignored" numbers*/
494             std::strcat(f1, "%lf"); /*and f1 to read the actual number in this pass of the loop*/
495             if (sscanf(str, f1, &d) != 1)
496             {
497                 gmx_fatal(FARGS, "Not enough elements for %s parameter (I need %d)", par, nvecs);
498             }
499             vec_params[i] = d;
500             std::strcat(f0, "%*s");
501         }
502     }
503     return vec_params;
504 }
505
506
507 static void init_edx(struct edix *edx)
508 {
509     edx->nr = 0;
510     snew(edx->x, 1);
511     snew(edx->anrs, 1);
512 }
513
514 static void filter2edx(struct edix *edx, int nindex, int index[], int ngro,
515                        int igro[], const rvec *x, const char* structure)
516 {
517 /* filter2edx copies coordinates from x to edx which are given in index
518  */
519
520     int pos, i;
521     int ix = edx->nr;
522     edx->nr += nindex;
523     srenew(edx->x, edx->nr);
524     srenew(edx->anrs, edx->nr);
525     for (i = 0; i < nindex; i++, ix++)
526     {
527         for (pos = 0; pos < ngro-1 && igro[pos] != index[i]; ++pos)
528         {
529         }
530         ;                                                            /*search element in igro*/
531         if (igro[pos] != index[i])
532         {
533             gmx_fatal(FARGS, "Couldn't find atom with index %d in structure %s", index[i], structure);
534         }
535         edx->anrs[ix] = index[i];
536         copy_rvec(x[pos], edx->x[ix]);
537     }
538 }
539
540 static void get_structure(const t_atoms *atoms, const char *IndexFile,
541                           const char *StructureFile, struct edix *edx, int nfit,
542                           int ifit[], int nav, int index[])
543 {
544     int     *igro; /*index corresponding to target or origin structure*/
545     int      ngro;
546     int      ntar;
547     rvec    *xtar;
548     char   * grpname;
549
550
551     ntar = read_conffile(StructureFile, &xtar);
552     printf("Select an index group of %d elements that corresponds to the atoms in the structure file %s\n",
553            ntar, StructureFile);
554     get_index(atoms, IndexFile, 1, &ngro, &igro, &grpname);
555     if (ngro != ntar)
556     {
557         gmx_fatal(FARGS, "You selected an index group with %d elements instead of %d", ngro, ntar);
558     }
559     init_edx(edx);
560     filter2edx(edx, nfit, ifit, ngro, igro, xtar, StructureFile);
561
562     /* If average and reference/fitting structure differ, append the average structure as well */
563     if (ifit != index) /*if fit structure is different append these coordinates, too -- don't mind duplicates*/
564     {
565         filter2edx(edx, nav, index, ngro, igro, xtar, StructureFile);
566     }
567 }
568
569 int gmx_make_edi(int argc, char *argv[])
570 {
571
572     static const char *desc[] = {
573         "[THISMODULE] generates an essential dynamics (ED) sampling input file to be used with [TT]mdrun[tt]",
574         "based on eigenvectors of a covariance matrix ([gmx-covar]) or from a",
575         "normal modes analysis ([gmx-nmeig]).",
576         "ED sampling can be used to manipulate the position along collective coordinates",
577         "(eigenvectors) of (biological) macromolecules during a simulation. Particularly,",
578         "it may be used to enhance the sampling efficiency of MD simulations by stimulating",
579         "the system to explore new regions along these collective coordinates. A number",
580         "of different algorithms are implemented to drive the system along the eigenvectors",
581         "([TT]-linfix[tt], [TT]-linacc[tt], [TT]-radfix[tt], [TT]-radacc[tt], [TT]-radcon[tt]),",
582         "to keep the position along a certain (set of) coordinate(s) fixed ([TT]-linfix[tt]),",
583         "or to only monitor the projections of the positions onto",
584         "these coordinates ([TT]-mon[tt]).[PAR]",
585         "References:[PAR]",
586         "A. Amadei, A.B.M. Linssen, B.L. de Groot, D.M.F. van Aalten and ",
587         "H.J.C. Berendsen; An efficient method for sampling the essential subspace ",
588         "of proteins., J. Biomol. Struct. Dyn. 13:615-626 (1996)[PAR]",
589         "B.L. de Groot, A. Amadei, D.M.F. van Aalten and H.J.C. Berendsen; ",
590         "Towards an exhaustive sampling of the configurational spaces of the ",
591         "two forms of the peptide hormone guanylin,",
592         "J. Biomol. Struct. Dyn. 13 : 741-751 (1996)[PAR]",
593         "B.L. de Groot, A.Amadei, R.M. Scheek, N.A.J. van Nuland and H.J.C. Berendsen; ",
594         "An extended sampling of the configurational space of HPr from E. coli",
595         "Proteins: Struct. Funct. Gen. 26: 314-322 (1996)",
596         "[PAR]You will be prompted for one or more index groups that correspond to the eigenvectors,",
597         "reference structure, target positions, etc.[PAR]",
598
599         "[TT]-mon[tt]: monitor projections of the coordinates onto selected eigenvectors.[PAR]",
600         "[TT]-linfix[tt]: perform fixed-step linear expansion along selected eigenvectors.[PAR]",
601         "[TT]-linacc[tt]: perform acceptance linear expansion along selected eigenvectors.",
602         "(steps in the desired directions will be accepted, others will be rejected).[PAR]",
603         "[TT]-radfix[tt]: perform fixed-step radius expansion along selected eigenvectors.[PAR]",
604         "[TT]-radacc[tt]: perform acceptance radius expansion along selected eigenvectors.",
605         "(steps in the desired direction will be accepted, others will be rejected).",
606         "[BB]Note:[bb] by default the starting MD structure will be taken as origin of the first",
607         "expansion cycle for radius expansion. If [TT]-ori[tt] is specified, you will be able",
608         "to read in a structure file that defines an external origin.[PAR]",
609         "[TT]-radcon[tt]: perform acceptance radius contraction along selected eigenvectors",
610         "towards a target structure specified with [TT]-tar[tt].[PAR]",
611         "NOTE: each eigenvector can be selected only once. [PAR]",
612         "[TT]-outfrq[tt]: frequency (in steps) of writing out projections etc. to [REF].xvg[ref] file[PAR]",
613         "[TT]-slope[tt]: minimal slope in acceptance radius expansion. A new expansion",
614         "cycle will be started if the spontaneous increase of the radius (in nm/step)",
615         "is less than the value specified.[PAR]",
616         "[TT]-maxedsteps[tt]: maximum number of steps per cycle in radius expansion",
617         "before a new cycle is started.[PAR]",
618         "Note on the parallel implementation: since ED sampling is a 'global' thing",
619         "(collective coordinates etc.), at least on the 'protein' side, ED sampling",
620         "is not very parallel-friendly from an implementation point of view. Because",
621         "parallel ED requires some extra communication, expect the performance to be",
622         "lower as in a free MD simulation, especially on a large number of ranks and/or",
623         "when the ED group contains a lot of atoms. [PAR]",
624         "Please also note that if your ED group contains more than a single protein,",
625         "then the [REF].tpr[ref] file must contain the correct PBC representation of the ED group.",
626         "Take a look on the initial RMSD from the reference structure, which is printed",
627         "out at the start of the simulation; if this is much higher than expected, one",
628         "of the ED molecules might be shifted by a box vector. [PAR]",
629         "All ED-related output of [TT]mdrun[tt] (specify with [TT]-eo[tt]) is written to a [REF].xvg[ref] file",
630         "as a function of time in intervals of OUTFRQ steps.[PAR]",
631         "[BB]Note[bb] that you can impose multiple ED constraints and flooding potentials in",
632         "a single simulation (on different molecules) if several [REF].edi[ref] files were concatenated",
633         "first. The constraints are applied in the order they appear in the [REF].edi[ref] file. ",
634         "Depending on what was specified in the [REF].edi[ref] input file, the output file contains for each ED dataset",
635         "",
636         " * the RMSD of the fitted molecule to the reference structure (for atoms involved in fitting prior to calculating the ED constraints)",
637         " * projections of the positions onto selected eigenvectors",
638         "",
639         "FLOODING:[PAR]",
640         "with [TT]-flood[tt], you can specify which eigenvectors are used to compute a flooding potential,",
641         "which will lead to extra forces expelling the structure out of the region described",
642         "by the covariance matrix. If you switch -restrain the potential is inverted and the structure",
643         "is kept in that region.",
644         "[PAR]",
645         "The origin is normally the average structure stored in the [TT]eigvec.trr[tt] file.",
646         "It can be changed with [TT]-ori[tt] to an arbitrary position in configuration space.",
647         "With [TT]-tau[tt], [TT]-deltaF0[tt], and [TT]-Eflnull[tt] you control the flooding behaviour.",
648         "Efl is the flooding strength, it is updated according to the rule of adaptive flooding.",
649         "Tau is the time constant of adaptive flooding, high [GRK]tau[grk] means slow adaption (i.e. growth). ",
650         "DeltaF0 is the flooding strength you want to reach after tau ps of simulation.",
651         "To use constant Efl set [TT]-tau[tt] to zero.",
652         "[PAR]",
653         "[TT]-alpha[tt] is a fudge parameter to control the width of the flooding potential. A value of 2 has been found",
654         "to give good results for most standard cases in flooding of proteins.",
655         "[GRK]alpha[grk] basically accounts for incomplete sampling, if you sampled further the width of the ensemble would",
656         "increase, this is mimicked by [GRK]alpha[grk] > 1.",
657         "For restraining, [GRK]alpha[grk] < 1 can give you smaller width in the restraining potential.",
658         "[PAR]",
659         "RESTART and FLOODING:",
660         "If you want to restart a crashed flooding simulation please find the values deltaF and Efl in",
661         "the output file and manually put them into the [REF].edi[ref] file under DELTA_F0 and EFL_NULL."
662     };
663
664     /* Save all the params in this struct and then save it in an edi file.
665      * ignoring fields nmass,massnrs,mass,tmass,nfit,fitnrs,edo
666      */
667     static t_edipar edi_params;
668
669     enum  {
670         evStepNr = evRADFIX + 1
671     };
672     static const char* evSelections[evNr]      = {nullptr, nullptr, nullptr, nullptr, nullptr, nullptr};
673     static const char* evOptions[evNr]         = {"-linfix", "-linacc", "-flood", "-radfix", "-radacc", "-radcon", "-mon"};
674     static const char* evParams[evStepNr]      = {nullptr, nullptr};
675     static const char* evStepOptions[evStepNr] = {"-linstep", "-accdir", "-not_used", "-radstep"};
676     static const char* ConstForceStr;
677     static real      * evStepList[evStepNr];
678     static real        radstep  = 0.0;
679     static real        deltaF0  = 150;
680     static real        deltaF   = 0;
681     static real        tau      = .1;
682     static real        constEfl = 0.0;
683     static real        alpha    = 1;
684     static int         eqSteps  = 0;
685     static int       * listen[evNr];
686     static real        T         = 300.0;
687     const real         kB        = 2.5 / 300.0; /* k_boltzmann in MD units */
688     static gmx_bool    bRestrain = FALSE;
689     static gmx_bool    bHesse    = FALSE;
690     static gmx_bool    bHarmonic = FALSE;
691     t_pargs            pa[]      = {
692         { "-mon", FALSE, etSTR, {&evSelections[evMON]},
693           "Indices of eigenvectors for projections of x (e.g. 1,2-5,9) or 1-100:10 means 1 11 21 31 ... 91" },
694         { "-linfix", FALSE, etSTR, {&evSelections[0]},
695           "Indices of eigenvectors for fixed increment linear sampling" },
696         { "-linacc", FALSE, etSTR, {&evSelections[1]},
697           "Indices of eigenvectors for acceptance linear sampling" },
698         { "-radfix", FALSE, etSTR, {&evSelections[3]},
699           "Indices of eigenvectors for fixed increment radius expansion" },
700         { "-radacc", FALSE, etSTR, {&evSelections[4]},
701           "Indices of eigenvectors for acceptance radius expansion" },
702         { "-radcon", FALSE, etSTR, {&evSelections[5]},
703           "Indices of eigenvectors for acceptance radius contraction" },
704         { "-flood",  FALSE, etSTR, {&evSelections[2]},
705           "Indices of eigenvectors for flooding"},
706         { "-outfrq", FALSE, etINT, {&edi_params.outfrq},
707           "Frequency (in steps) of writing output in [REF].xvg[ref] file" },
708         { "-slope", FALSE, etREAL, { &edi_params.slope},
709           "Minimal slope in acceptance radius expansion"},
710         { "-linstep", FALSE, etSTR, {&evParams[0]},
711           "Stepsizes (nm/step) for fixed increment linear sampling (put in quotes! \"1.0 2.3 5.1 -3.1\")"},
712         { "-accdir", FALSE, etSTR, {&evParams[1]},
713           "Directions for acceptance linear sampling - only sign counts! (put in quotes! \"-1 +1 -1.1\")"},
714         { "-radstep", FALSE, etREAL, {&radstep},
715           "Stepsize (nm/step) for fixed increment radius expansion"},
716         { "-maxedsteps", FALSE, etINT, {&edi_params.maxedsteps},
717           "Maximum number of steps per cycle" },
718         { "-eqsteps", FALSE, etINT, {&eqSteps},
719           "Number of steps to run without any perturbations "},
720         { "-deltaF0", FALSE, etREAL, {&deltaF0},
721           "Target destabilization energy for flooding"},
722         { "-deltaF", FALSE, etREAL, {&deltaF},
723           "Start deltaF with this parameter - default 0, nonzero values only needed for restart"},
724         { "-tau", FALSE, etREAL, {&tau},
725           "Coupling constant for adaption of flooding strength according to deltaF0, 0 = infinity i.e. constant flooding strength"},
726         { "-Eflnull", FALSE, etREAL, {&constEfl},
727           "The starting value of the flooding strength. The flooding strength is updated "
728           "according to the adaptive flooding scheme. For a constant flooding strength use [TT]-tau[tt] 0. "},
729         { "-T", FALSE, etREAL, {&T},
730           "T is temperature, the value is needed if you want to do flooding "},
731         { "-alpha", FALSE, etREAL, {&alpha},
732           "Scale width of gaussian flooding potential with alpha^2 "},
733         { "-restrain", FALSE, etBOOL, {&bRestrain},
734           "Use the flooding potential with inverted sign -> effects as quasiharmonic restraining potential"},
735         { "-hessian", FALSE, etBOOL, {&bHesse},
736           "The eigenvectors and eigenvalues are from a Hessian matrix"},
737         { "-harmonic", FALSE, etBOOL, {&bHarmonic},
738           "The eigenvalues are interpreted as spring constant"},
739         { "-constF", FALSE, etSTR, {&ConstForceStr},
740           "Constant force flooding: manually set the forces for the eigenvectors selected with -flood "
741           "(put in quotes! \"1.0 2.3 5.1 -3.1\"). No other flooding parameters are needed when specifying the forces directly."}
742     };
743 #define NPA asize(pa)
744
745     rvec             *xref1;
746     int               nvec1, *eignr1 = nullptr;
747     rvec             *xav1, **eigvec1 = nullptr;
748     t_atoms          *atoms = nullptr;
749     int               nav; /* Number of atoms in the average structure */
750     char             *grpname;
751     const char       *indexfile;
752     int               i;
753     int              *index, *ifit;
754     int               nfit;              /* Number of atoms in the reference/fit structure */
755     int               ev_class;          /* parameter _class i.e. evMON, evRADFIX etc. */
756     int               nvecs;
757     real             *eigval1 = nullptr; /* in V3.3 this is parameter of read_eigenvectors */
758
759     const char       *EdiFile;
760     const char       *TargetFile;
761     const char       *OriginFile;
762     const char       *EigvecFile;
763
764     gmx_output_env_t *oenv;
765
766     /*to read topology file*/
767     t_topology  top;
768     int         ePBC;
769     matrix      topbox;
770     rvec       *xtop;
771     gmx_bool    bFit1;
772
773     t_filenm    fnm[] = {
774         { efTRN, "-f",    "eigenvec",    ffREAD  },
775         { efXVG, "-eig",  "eigenval",    ffOPTRD },
776         { efTPS, nullptr,    nullptr,          ffREAD },
777         { efNDX, nullptr,    nullptr,  ffOPTRD },
778         { efSTX, "-tar", "target", ffOPTRD},
779         { efSTX, "-ori", "origin", ffOPTRD},
780         { efEDI, "-o", "sam", ffWRITE }
781     };
782 #define NFILE asize(fnm)
783     edi_params.outfrq = 100; edi_params.slope = 0.0; edi_params.maxedsteps = 0;
784     if (!parse_common_args(&argc, argv, 0,
785                            NFILE, fnm, NPA, pa, asize(desc), desc, 0, nullptr, &oenv))
786     {
787         return 0;
788     }
789
790     indexfile       = ftp2fn_null(efNDX, NFILE, fnm);
791     EdiFile         = ftp2fn(efEDI, NFILE, fnm);
792     TargetFile      = opt2fn_null("-tar", NFILE, fnm);
793     OriginFile      = opt2fn_null("-ori", NFILE, fnm);
794
795
796     for (ev_class = 0; ev_class < evNr; ++ev_class)
797     {
798         if (opt2parg_bSet(evOptions[ev_class], NPA, pa))
799         {
800             /*get list of eigenvectors*/
801             nvecs = sscan_list(&(listen[ev_class]), opt2parg_str(evOptions[ev_class], NPA, pa), evOptions[ev_class]);
802             if (ev_class < evStepNr-2)
803             {
804                 /*if apropriate get list of stepsizes for these eigenvectors*/
805                 if (opt2parg_bSet(evStepOptions[ev_class], NPA, pa))
806                 {
807                     evStepList[ev_class] =
808                         scan_vecparams(opt2parg_str(evStepOptions[ev_class], NPA, pa), evStepOptions[ev_class], nvecs);
809                 }
810                 else   /*if list is not given fill with zeros */
811                 {
812                     snew(evStepList[ev_class], nvecs);
813                     for (i = 0; i < nvecs; i++)
814                     {
815                         evStepList[ev_class][i] = 0.0;
816                     }
817                 }
818             }
819             else if (ev_class == evRADFIX)
820             {
821                 snew(evStepList[ev_class], nvecs);
822                 for (i = 0; i < nvecs; i++)
823                 {
824                     evStepList[ev_class][i] = radstep;
825                 }
826             }
827             else if (ev_class == evFLOOD)
828             {
829                 snew(evStepList[ev_class], nvecs);
830
831                 /* Are we doing constant force flooding? In that case, we read in
832                  * the fproj values from the command line */
833                 if (opt2parg_bSet("-constF", NPA, pa))
834                 {
835                     evStepList[ev_class] = scan_vecparams(opt2parg_str("-constF", NPA, pa), "-constF", nvecs);
836                 }
837             }
838             else
839             {
840             };   /*to avoid ambiguity   */
841         }
842         else     /* if there are no eigenvectors for this option set list to zero */
843         {
844             listen[ev_class] = nullptr;
845             snew(listen[ev_class], 1);
846             listen[ev_class][0] = 0;
847         }
848     }
849
850     /* print the interpreted list of eigenvectors - to give some feedback*/
851     for (ev_class = 0; ev_class < evNr; ++ev_class)
852     {
853         printf("Eigenvector list %7s consists of the indices: ", evOptions[ev_class]);
854         i = 0;
855         while (listen[ev_class][i])
856         {
857             printf("%d ", listen[ev_class][i++]);
858         }
859         printf("\n");
860     }
861
862     EigvecFile = opt2fn("-f", NFILE, fnm);
863
864     /*read eigenvectors from eigvec.trr*/
865     read_eigenvectors(EigvecFile, &nav, &bFit1,
866                       &xref1, &edi_params.fitmas, &xav1, &edi_params.pcamas, &nvec1, &eignr1, &eigvec1, &eigval1);
867
868     read_tps_conf(ftp2fn(efTPS, NFILE, fnm),
869                   &top, &ePBC, &xtop, nullptr, topbox, 0);
870     atoms = &top.atoms;
871
872
873     printf("\nSelect an index group of %d elements that corresponds to the eigenvectors\n", nav);
874     get_index(atoms, indexfile, 1, &i, &index, &grpname); /*if indexfile != NULL parameter 'atoms' is ignored */
875     if (i != nav)
876     {
877         gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
878                   i, nav);
879     }
880     printf("\n");
881
882
883     if (xref1 == nullptr)
884     {
885         if (bFit1)
886         {
887             /* if g_covar used different coordinate groups to fit and to do the PCA */
888             printf("\nNote: the structure in %s should be the same\n"
889                    "      as the one used for the fit in g_covar\n", ftp2fn(efTPS, NFILE, fnm));
890             printf("\nSelect the index group that was used for the least squares fit in g_covar\n");
891         }
892         else
893         {
894             printf("\nNote: Apparently no fitting was done in g_covar.\n"
895                    "      However, you need to select a reference group for fitting in mdrun\n");
896         }
897         get_index(atoms, indexfile, 1, &nfit, &ifit, &grpname);
898         snew(xref1, nfit);
899         for (i = 0; i < nfit; i++)
900         {
901             copy_rvec(xtop[ifit[i]], xref1[i]);
902         }
903     }
904     else
905     {
906         nfit = nav;
907         ifit = index;
908     }
909
910     if (opt2parg_bSet("-constF", NPA, pa))
911     {
912         /* Constant force flooding is special: Most of the normal flooding
913          * options are not needed. */
914         edi_params.flood.bConstForce = TRUE;
915     }
916     else
917     {
918         /* For normal flooding read eigenvalues and store them in evSteplist[evFLOOD] */
919
920         if (listen[evFLOOD][0] != 0)
921         {
922             read_eigenvalues(listen[evFLOOD], opt2fn("-eig", NFILE, fnm), evStepList[evFLOOD], bHesse, kB*T, nav);
923         }
924
925         edi_params.flood.tau       = tau;
926         edi_params.flood.deltaF0   = deltaF0;
927         edi_params.flood.deltaF    = deltaF;
928         edi_params.presteps        = eqSteps;
929         edi_params.flood.kT        = kB*T;
930         edi_params.flood.bHarmonic = bHarmonic;
931         if (bRestrain)
932         {
933             /* Trick: invert sign of Efl and alpha2 then this will give the same sign in the exponential and inverted sign outside */
934             edi_params.flood.constEfl = -constEfl;
935             edi_params.flood.alpha2   = -gmx::square(alpha);
936         }
937         else
938         {
939             edi_params.flood.constEfl = constEfl;
940             edi_params.flood.alpha2   = gmx::square(alpha);
941         }
942     }
943
944     edi_params.ned = nav;
945
946     /*number of system atoms  */
947     edi_params.nini = atoms->nr;
948
949
950     /*store reference and average structure in edi_params*/
951     make_t_edx(&edi_params.sref, nfit, xref1, ifit );
952     make_t_edx(&edi_params.sav, nav, xav1, index);
953
954
955     /* Store target positions in edi_params */
956     if (opt2bSet("-tar", NFILE, fnm))
957     {
958         if (0 != listen[evFLOOD][0])
959         {
960             fprintf(stderr, "\nNote: Providing a TARGET structure has no effect when using flooding.\n"
961                     "      You may want to use -ori to define the flooding potential center.\n\n");
962         }
963         get_structure(atoms, indexfile, TargetFile, &edi_params.star, nfit, ifit, nav, index);
964     }
965     else
966     {
967         make_t_edx(&edi_params.star, 0, nullptr, index);
968     }
969
970     /* Store origin positions */
971     if (opt2bSet("-ori", NFILE, fnm))
972     {
973         get_structure(atoms, indexfile, OriginFile, &edi_params.sori, nfit, ifit, nav, index);
974     }
975     else
976     {
977         make_t_edx(&edi_params.sori, 0, nullptr, index);
978     }
979
980     /* Write edi-file */
981     write_the_whole_thing(gmx_ffopen(EdiFile, "w"), &edi_params, eigvec1, nvec1, listen, evStepList);
982
983     return 0;
984 }