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