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