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