2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012, 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.
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.
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.
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.
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.
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.
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.
53 #include "gmx_fatal.h"
56 #include "gromacs/fileio/futil.h"
58 #include "gromacs/fileio/pdbio.h"
59 #include "gromacs/fileio/confio.h"
60 #include "gromacs/fileio/tpxio.h"
61 #include "gromacs/fileio/matio.h"
75 gmx_bool bConstForce; /* Do constant force flooding instead of
76 evaluating a flooding potential */
85 /* This type is for the average, reference, target, and origin structure */
88 int nr; /* number of atoms this structure contains */
89 int *anrs; /* atom index numbers */
90 rvec *x; /* positions */
91 real *sqrtm; /* sqrt of the masses used for mass-
92 * weighting of analysis */
98 int nini; /* total Nr of atoms */
99 gmx_bool fitmas; /* true if trans fit with cm */
100 gmx_bool pcamas; /* true if mass-weighted PCA */
101 int presteps; /* number of steps to run without any
102 * perturbations ... just monitoring */
103 int outfrq; /* freq (in steps) of writing to edo */
104 int maxedsteps; /* max nr of steps per cycle */
105 struct edix sref; /* reference positions, to these fitting
107 struct edix sav; /* average positions */
108 struct edix star; /* target positions */
109 struct edix sori; /* origin positions */
110 real slope; /* minimal slope in acceptance radexp */
111 int ned; /* Nr of atoms in essdyn buffer */
112 t_edflood flood; /* parameters especially for flooding */
117 void make_t_edx(struct edix *edx, int natoms, rvec *pos, atom_id index[])
124 void write_t_edx(FILE *fp, struct edix edx, const char *comment)
126 /*here we copy only the pointers into the t_edx struct
127 no data is copied and edx.box is ignored */
129 fprintf(fp, "#%s \n %d \n", comment, edx.nr);
130 for (i = 0; i < edx.nr; i++)
132 fprintf(fp, "%d %f %f %f\n", (edx.anrs)[i]+1, (edx.x)[i][XX], (edx.x)[i][YY], (edx.x)[i][ZZ]);
136 int sscan_list(int *list[], const char *str, const char *listname)
138 /*this routine scans a string of the form 1,3-6,9 and returns the
139 selected numbers (in this case 1 3 4 5 6 9) in NULL-terminated array of integers.
140 memory for this list will be allocated in this routine -- sscan_list expects *list to
143 listname is a string used in the errormessage*/
148 char *pos, *startpos, *step;
151 /*enums to define the different lexical stati */
153 sBefore, sNumber, sMinus, sRange, sZero, sSmaller, sError, sSteppedRange
156 int status = sBefore; /*status of the deterministic automat to scan str */
160 char *start = NULL; /*holds the string of the number behind a ','*/
161 char *end = NULL; /*holds the string of the number behind a '-' */
163 int nvecs = 0; /* counts the number of vectors in the list*/
175 while ((c = *pos) != 0)
179 /* expect a number */
180 case sBefore: if (isdigit(c))
191 /* have read a number, expect ',' or '-' */
192 case sNumber: if (c == ',')
195 srenew(*list, nvecs+1);
196 (*list)[nvecs++] = number = strtol(start, NULL, 10);
206 status = sMinus; break;
217 /* have read a '-' -> expect a number */
222 status = sRange; break;
234 status = sError; break;
248 /* have read the number after a minus, expect ',' or ':' */
253 end_number = strtol(end, NULL, 10);
254 number = strtol(start, NULL, 10);
258 status = sZero; break;
260 if (end_number <= number)
262 status = sSmaller; break;
264 srenew(*list, nvecs+end_number-number+1);
267 istep = strtol(step, NULL, 10);
274 for (i = number; i <= end_number; i += istep)
276 (*list)[nvecs++] = i;
282 status = sSteppedRange;
294 /* format error occured */
296 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %d with char %c", listname, pos-startpos, *(pos-1));
298 /* logical error occured */
300 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %d: eigenvector 0 is not valid", listname, pos-startpos);
303 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 ++pos; /* read next character */
307 } /*scanner has finished */
309 /* append zero to list of eigenvectors */
310 srenew(*list, nvecs+1);
316 void write_eigvec(FILE* fp, int natoms, int eig_list[], rvec** eigvecs, int nvec, const char *grouptitle, real steps[])
318 /* eig_list is a zero-terminated list of indices into the eigvecs array.
319 eigvecs are coordinates of eigenvectors
320 grouptitle to write in the comment line
321 steps -- array with stepsizes for evLINFIX, evLINACC and evRADACC
324 int n = 0, i; rvec x;
326 while (eig_list[n++])
328 ; /*count selected eigenvecs*/
331 fprintf(fp, "# NUMBER OF EIGENVECTORS + %s\n %d\n", grouptitle, n-1);
333 /* write list of eigenvector indicess */
334 for (n = 0; eig_list[n]; n++)
338 fprintf(fp, "%8d %g\n", eig_list[n], steps[n]);
342 fprintf(fp, "%8d %g\n", eig_list[n], 1.0);
347 /* dump coordinates of the selected eigenvectors */
351 for (i = 0; i < natoms; i++)
353 if (eig_list[n] > nvec)
355 gmx_fatal(FARGS, "Selected eigenvector %d is higher than maximum number %d of available eigenvectors", eig_list[n], nvec);
357 copy_rvec(eigvecs[eig_list[n]-1][i], x);
359 fprintf(fp, "%8.5f %8.5f %8.5f\n", x[XX], x[YY], x[ZZ]);
366 /*enum referring to the different lists of eigenvectors*/
368 evLINFIX, evLINACC, evFLOOD, evRADFIX, evRADACC, evRADCON, evMON, evNr
374 void write_the_whole_thing(FILE* fp, t_edipar *edpars, rvec** eigvecs,
375 int nvec, int *eig_listen[], real* evStepList[])
380 fprintf(fp, "#MAGIC\n %d \n#NINI\n %d\n#FITMAS\n %d\n#ANALYSIS_MAS\n %d\n",
381 MAGIC, edpars->nini, edpars->fitmas, edpars->pcamas);
382 fprintf(fp, "#OUTFRQ\n %d\n#MAXLEN\n %d\n#SLOPECRIT\n %f\n",
383 edpars->outfrq, edpars->maxedsteps, edpars->slope);
384 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",
385 edpars->presteps, edpars->flood.deltaF0, edpars->flood.deltaF, edpars->flood.tau, edpars->flood.constEfl,
386 edpars->flood.alpha2, edpars->flood.kT, edpars->flood.bHarmonic, edpars->flood.bConstForce);
388 /* Average and reference positions */
389 write_t_edx(fp, edpars->sref, "NREF, XREF");
390 write_t_edx(fp, edpars->sav, "NAV, XAV");
394 write_eigvec(fp, edpars->ned, eig_listen[evMON], eigvecs, nvec, "COMPONENTS GROUP 1", NULL);
395 write_eigvec(fp, edpars->ned, eig_listen[evLINFIX], eigvecs, nvec, "COMPONENTS GROUP 2", evStepList[evLINFIX]);
396 write_eigvec(fp, edpars->ned, eig_listen[evLINACC], eigvecs, nvec, "COMPONENTS GROUP 3", evStepList[evLINACC]);
397 write_eigvec(fp, edpars->ned, eig_listen[evRADFIX], eigvecs, nvec, "COMPONENTS GROUP 4", evStepList[evRADFIX]);
398 write_eigvec(fp, edpars->ned, eig_listen[evRADACC], eigvecs, nvec, "COMPONENTS GROUP 5", NULL);
399 write_eigvec(fp, edpars->ned, eig_listen[evRADCON], eigvecs, nvec, "COMPONENTS GROUP 6", NULL);
400 write_eigvec(fp, edpars->ned, eig_listen[evFLOOD], eigvecs, nvec, "COMPONENTS GROUP 7", evStepList[evFLOOD]);
403 /*Target and Origin positions */
404 write_t_edx(fp, edpars->star, "NTARGET, XTARGET");
405 write_t_edx(fp, edpars->sori, "NORIGIN, XORIGIN");
408 int read_conffile(const char *confin, char *title, rvec *x[])
410 /* read coordinates out of STX file */
414 printf("read coordnumber from file %s\n", confin);
415 get_stx_coordnum(confin, &natoms);
416 printf("number of coordinates in file %d\n", natoms);
417 /* if (natoms != ncoords)
418 gmx_fatal(FARGS,"number of coordinates in coordinate file (%s, %d)\n"
419 " does not match topology (= %d)",
420 confin,natoms,ncoords);
422 /* make space for coordinates and velocities */
423 init_t_atoms(&confat, natoms, FALSE);
425 read_stx_conf(confin, title, &confat, *x, NULL, NULL, box);
430 void read_eigenvalues(int vecs[], const char *eigfile, real values[],
431 gmx_bool bHesse, real kT)
436 neig = read_xvg(eigfile, &eigval, &nrow);
438 fprintf(stderr, "Read %d eigenvalues\n", neig);
439 for (i = bHesse ? 6 : 0; i < neig; i++)
441 if (eigval[1][i] < -0.001 && bHesse)
444 "WARNING: The Hessian Matrix has negative eigenvalue %f, we set it to zero (no flooding in this direction)\n\n", eigval[1][i]);
447 if (eigval[1][i] < 0)
454 for (i = 0; vecs[i]; i++)
458 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");
460 values[i] = eigval[1][vecs[i]-1]/kT;
465 for (i = 0; vecs[i]; i++)
467 if (vecs[i] > (neig-6))
469 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");
471 values[i] = 1/eigval[1][vecs[i]-1];
475 for (i = 0; i < nrow; i++)
483 static real *scan_vecparams(const char *str, const char * par, int nvecs)
485 char f0[256], f1[256]; /*format strings adapted every pass of the loop*/
490 snew(vec_params, nvecs);
494 for (i = 0; (i < nvecs); i++)
496 strcpy(f1, f0); /*f0 is the format string for the "to-be-ignored" numbers*/
497 strcat(f1, "%lf"); /*and f1 to read the actual number in this pass of the loop*/
498 if (sscanf(str, f1, &d) != 1)
500 gmx_fatal(FARGS, "Not enough elements for %s parameter (I need %d)", par, nvecs);
511 void init_edx(struct edix *edx)
518 void filter2edx(struct edix *edx, int nindex, atom_id index[], int ngro,
519 atom_id igro[], rvec *x, const char* structure)
521 /* filter2edx copies coordinates from x to edx which are given in index
527 srenew(edx->x, edx->nr);
528 srenew(edx->anrs, edx->nr);
529 for (i = 0; i < nindex; i++, ix++)
531 for (pos = 0; pos < ngro-1 && igro[pos] != index[i]; ++pos)
534 ; /*search element in igro*/
535 if (igro[pos] != index[i])
537 gmx_fatal(FARGS, "Couldn't find atom with index %d in structure %s", index[i], structure);
539 edx->anrs[ix] = index[i];
540 copy_rvec(x[pos], edx->x[ix]);
544 void get_structure(t_atoms *atoms, const char *IndexFile,
545 const char *StructureFile, struct edix *edx, int nfit,
546 atom_id ifit[], int nav, atom_id index[])
548 atom_id *igro; /*index corresponding to target or origin structure*/
556 ntar = read_conffile(StructureFile, title, &xtar);
557 printf("Select an index group of %d elements that corresponds to the atoms in the structure file %s\n",
558 ntar, StructureFile);
559 get_index(atoms, IndexFile, 1, &ngro, &igro, &grpname);
562 gmx_fatal(FARGS, "You selected an index group with %d elements instead of %d", ngro, ntar);
565 filter2edx(edx, nfit, ifit, ngro, igro, xtar, StructureFile);
567 /* If average and reference/fitting structure differ, append the average structure as well */
568 if (ifit != index) /*if fit structure is different append these coordinates, too -- don't mind duplicates*/
570 filter2edx(edx, nav, index, ngro, igro, xtar, StructureFile);
574 int gmx_make_edi(int argc, char *argv[])
577 static const char *desc[] = {
578 "[TT]make_edi[tt] generates an essential dynamics (ED) sampling input file to be used with [TT]mdrun[tt]",
579 "based on eigenvectors of a covariance matrix ([TT]g_covar[tt]) or from a",
580 "normal modes analysis ([TT]g_nmeig[tt]).",
581 "ED sampling can be used to manipulate the position along collective coordinates",
582 "(eigenvectors) of (biological) macromolecules during a simulation. Particularly,",
583 "it may be used to enhance the sampling efficiency of MD simulations by stimulating",
584 "the system to explore new regions along these collective coordinates. A number",
585 "of different algorithms are implemented to drive the system along the eigenvectors",
586 "([TT]-linfix[tt], [TT]-linacc[tt], [TT]-radfix[tt], [TT]-radacc[tt], [TT]-radcon[tt]),",
587 "to keep the position along a certain (set of) coordinate(s) fixed ([TT]-linfix[tt]),",
588 "or to only monitor the projections of the positions onto",
589 "these coordinates ([TT]-mon[tt]).[PAR]",
591 "A. Amadei, A.B.M. Linssen, B.L. de Groot, D.M.F. van Aalten and ",
592 "H.J.C. Berendsen; An efficient method for sampling the essential subspace ",
593 "of proteins., J. Biomol. Struct. Dyn. 13:615-626 (1996)[BR]",
594 "B.L. de Groot, A. Amadei, D.M.F. van Aalten and H.J.C. Berendsen; ",
595 "Towards an exhaustive sampling of the configurational spaces of the ",
596 "two forms of the peptide hormone guanylin,",
597 "J. Biomol. Struct. Dyn. 13 : 741-751 (1996)[BR]",
598 "B.L. de Groot, A.Amadei, R.M. Scheek, N.A.J. van Nuland and H.J.C. Berendsen; ",
599 "An extended sampling of the configurational space of HPr from E. coli",
600 "Proteins: Struct. Funct. Gen. 26: 314-322 (1996)",
601 "[PAR]You will be prompted for one or more index groups that correspond to the eigenvectors,",
602 "reference structure, target positions, etc.[PAR]",
604 "[TT]-mon[tt]: monitor projections of the coordinates onto selected eigenvectors.[PAR]",
605 "[TT]-linfix[tt]: perform fixed-step linear expansion along selected eigenvectors.[PAR]",
606 "[TT]-linacc[tt]: perform acceptance linear expansion along selected eigenvectors.",
607 "(steps in the desired directions will be accepted, others will be rejected).[PAR]",
608 "[TT]-radfix[tt]: perform fixed-step radius expansion along selected eigenvectors.[PAR]",
609 "[TT]-radacc[tt]: perform acceptance radius expansion along selected eigenvectors.",
610 "(steps in the desired direction will be accepted, others will be rejected).",
611 "[BB]Note:[bb] by default the starting MD structure will be taken as origin of the first",
612 "expansion cycle for radius expansion. If [TT]-ori[tt] is specified, you will be able",
613 "to read in a structure file that defines an external origin.[PAR]",
614 "[TT]-radcon[tt]: perform acceptance radius contraction along selected eigenvectors",
615 "towards a target structure specified with [TT]-tar[tt].[PAR]",
616 "NOTE: each eigenvector can be selected only once. [PAR]",
617 "[TT]-outfrq[tt]: frequency (in steps) of writing out projections etc. to [TT].xvg[tt] file[PAR]",
618 "[TT]-slope[tt]: minimal slope in acceptance radius expansion. A new expansion",
619 "cycle will be started if the spontaneous increase of the radius (in nm/step)",
620 "is less than the value specified.[PAR]",
621 "[TT]-maxedsteps[tt]: maximum number of steps per cycle in radius expansion",
622 "before a new cycle is started.[PAR]",
623 "Note on the parallel implementation: since ED sampling is a 'global' thing",
624 "(collective coordinates etc.), at least on the 'protein' side, ED sampling",
625 "is not very parallel-friendly from an implementation point of view. Because",
626 "parallel ED requires some extra communication, expect the performance to be",
627 "lower as in a free MD simulation, especially on a large number of nodes and/or",
628 "when the ED group contains a lot of atoms. [PAR]",
629 "Please also note that if your ED group contains more than a single protein,",
630 "then the [TT].tpr[tt] file must contain the correct PBC representation of the ED group.",
631 "Take a look on the initial RMSD from the reference structure, which is printed",
632 "out at the start of the simulation; if this is much higher than expected, one",
633 "of the ED molecules might be shifted by a box vector. [PAR]",
634 "All ED-related output of [TT]mdrun[tt] (specify with [TT]-eo[tt]) is written to a [TT].xvg[tt] file",
635 "as a function of time in intervals of OUTFRQ steps.[PAR]",
636 "[BB]Note[bb] that you can impose multiple ED constraints and flooding potentials in",
637 "a single simulation (on different molecules) if several [TT].edi[tt] files were concatenated",
638 "first. The constraints are applied in the order they appear in the [TT].edi[tt] file. ",
639 "Depending on what was specified in the [TT].edi[tt] input file, the output file contains for each ED dataset[PAR]",
640 "[TT]*[tt] the RMSD of the fitted molecule to the reference structure (for atoms involved in fitting prior to calculating the ED constraints)[BR]",
641 "[TT]*[tt] projections of the positions onto selected eigenvectors[BR]",
644 "with [TT]-flood[tt], you can specify which eigenvectors are used to compute a flooding potential,",
645 "which will lead to extra forces expelling the structure out of the region described",
646 "by the covariance matrix. If you switch -restrain the potential is inverted and the structure",
647 "is kept in that region.",
649 "The origin is normally the average structure stored in the [TT]eigvec.trr[tt] file.",
650 "It can be changed with [TT]-ori[tt] to an arbitrary position in configuration space.",
651 "With [TT]-tau[tt], [TT]-deltaF0[tt], and [TT]-Eflnull[tt] you control the flooding behaviour.",
652 "Efl is the flooding strength, it is updated according to the rule of adaptive flooding.",
653 "Tau is the time constant of adaptive flooding, high [GRK]tau[grk] means slow adaption (i.e. growth). ",
654 "DeltaF0 is the flooding strength you want to reach after tau ps of simulation.",
655 "To use constant Efl set [TT]-tau[tt] to zero.",
657 "[TT]-alpha[tt] is a fudge parameter to control the width of the flooding potential. A value of 2 has been found",
658 "to give good results for most standard cases in flooding of proteins.",
659 "[GRK]alpha[grk] basically accounts for incomplete sampling, if you sampled further the width of the ensemble would",
660 "increase, this is mimicked by [GRK]alpha[grk] > 1.",
661 "For restraining, [GRK]alpha[grk] < 1 can give you smaller width in the restraining potential.",
663 "RESTART and FLOODING:",
664 "If you want to restart a crashed flooding simulation please find the values deltaF and Efl in",
665 "the output file and manually put them into the [TT].edi[tt] file under DELTA_F0 and EFL_NULL."
668 /* Save all the params in this struct and then save it in an edi file.
669 * ignoring fields nmass,massnrs,mass,tmass,nfit,fitnrs,edo
671 static t_edipar edi_params;
674 evStepNr = evRADFIX + 1
676 static const char* evSelections[evNr] = {NULL, NULL, NULL, NULL, NULL, NULL};
677 static const char* evOptions[evNr] = {"-linfix", "-linacc", "-flood", "-radfix", "-radacc", "-radcon", "-mon"};
678 static const char* evParams[evStepNr] = {NULL, NULL};
679 static const char* evStepOptions[evStepNr] = {"-linstep", "-accdir", "-not_used", "-radstep"};
680 static const char* ConstForceStr;
681 static real * evStepList[evStepNr];
682 static real radfix = 0.0;
683 static real deltaF0 = 150;
684 static real deltaF = 0;
685 static real tau = .1;
686 static real constEfl = 0.0;
687 static real alpha = 1;
688 static int eqSteps = 0;
689 static int * listen[evNr];
690 static real T = 300.0;
691 const real kB = 2.5 / 300.0; /* k_boltzmann in MD units */
692 static gmx_bool bRestrain = FALSE;
693 static gmx_bool bHesse = FALSE;
694 static gmx_bool bHarmonic = FALSE;
696 { "-mon", FALSE, etSTR, {&evSelections[evMON]},
697 "Indices of eigenvectors for projections of x (e.g. 1,2-5,9) or 1-100:10 means 1 11 21 31 ... 91" },
698 { "-linfix", FALSE, etSTR, {&evSelections[0]},
699 "Indices of eigenvectors for fixed increment linear sampling" },
700 { "-linacc", FALSE, etSTR, {&evSelections[1]},
701 "Indices of eigenvectors for acceptance linear sampling" },
702 { "-radfix", FALSE, etSTR, {&evSelections[3]},
703 "Indices of eigenvectors for fixed increment radius expansion" },
704 { "-radacc", FALSE, etSTR, {&evSelections[4]},
705 "Indices of eigenvectors for acceptance radius expansion" },
706 { "-radcon", FALSE, etSTR, {&evSelections[5]},
707 "Indices of eigenvectors for acceptance radius contraction" },
708 { "-flood", FALSE, etSTR, {&evSelections[2]},
709 "Indices of eigenvectors for flooding"},
710 { "-outfrq", FALSE, etINT, {&edi_params.outfrq},
711 "Freqency (in steps) of writing output in [TT].xvg[tt] file" },
712 { "-slope", FALSE, etREAL, { &edi_params.slope},
713 "Minimal slope in acceptance radius expansion"},
714 { "-linstep", FALSE, etSTR, {&evParams[0]},
715 "Stepsizes (nm/step) for fixed increment linear sampling (put in quotes! \"1.0 2.3 5.1 -3.1\")"},
716 { "-accdir", FALSE, etSTR, {&evParams[1]},
717 "Directions for acceptance linear sampling - only sign counts! (put in quotes! \"-1 +1 -1.1\")"},
718 { "-radstep", FALSE, etREAL, {&radfix},
719 "Stepsize (nm/step) for fixed increment radius expansion"},
720 { "-maxedsteps", FALSE, etINT, {&edi_params.maxedsteps},
721 "Maximum number of steps per cycle" },
722 { "-eqsteps", FALSE, etINT, {&eqSteps},
723 "Number of steps to run without any perturbations "},
724 { "-deltaF0", FALSE, etREAL, {&deltaF0},
725 "Target destabilization energy for flooding"},
726 { "-deltaF", FALSE, etREAL, {&deltaF},
727 "Start deltaF with this parameter - default 0, nonzero values only needed for restart"},
728 { "-tau", FALSE, etREAL, {&tau},
729 "Coupling constant for adaption of flooding strength according to deltaF0, 0 = infinity i.e. constant flooding strength"},
730 { "-Eflnull", FALSE, etREAL, {&constEfl},
731 "The starting value of the flooding strength. The flooding strength is updated "
732 "according to the adaptive flooding scheme. For a constant flooding strength use [TT]-tau[tt] 0. "},
733 { "-T", FALSE, etREAL, {&T},
734 "T is temperature, the value is needed if you want to do flooding "},
735 { "-alpha", FALSE, etREAL, {&alpha},
736 "Scale width of gaussian flooding potential with alpha^2 "},
737 { "-restrain", FALSE, etBOOL, {&bRestrain},
738 "Use the flooding potential with inverted sign -> effects as quasiharmonic restraining potential"},
739 { "-hessian", FALSE, etBOOL, {&bHesse},
740 "The eigenvectors and eigenvalues are from a Hessian matrix"},
741 { "-harmonic", FALSE, etBOOL, {&bHarmonic},
742 "The eigenvalues are interpreted as spring constant"},
743 { "-constF", FALSE, etSTR, {&ConstForceStr},
744 "Constant force flooding: manually set the forces for the eigenvectors selected with -flood "
745 "(put in quotes! \"1.0 2.3 5.1 -3.1\"). No other flooding parameters are needed when specifying the forces directly."}
747 #define NPA asize(pa)
750 int nvec1, *eignr1 = NULL;
751 rvec *xav1, **eigvec1 = NULL;
752 t_atoms *atoms = NULL;
753 int nav; /* Number of atoms in the average structure */
755 const char *indexfile;
757 atom_id *index, *ifit;
758 int nfit; /* Number of atoms in the reference/fit structure */
759 int ev_class; /* parameter _class i.e. evMON, evRADFIX etc. */
761 real *eigval1 = NULL; /* in V3.3 this is parameter of read_eigenvectors */
764 const char *TargetFile;
765 const char *OriginFile;
766 const char *EigvecFile;
770 /*to read topology file*/
776 gmx_bool bTop, bFit1;
779 { efTRN, "-f", "eigenvec", ffREAD },
780 { efXVG, "-eig", "eigenval", ffOPTRD },
781 { efTPS, NULL, NULL, ffREAD },
782 { efNDX, NULL, NULL, ffOPTRD },
783 { efSTX, "-tar", "target", ffOPTRD},
784 { efSTX, "-ori", "origin", ffOPTRD},
785 { efEDI, "-o", "sam", ffWRITE }
787 #define NFILE asize(fnm)
788 edi_params.outfrq = 100; edi_params.slope = 0.0; edi_params.maxedsteps = 0;
789 if (!parse_common_args(&argc, argv, 0,
790 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
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);
801 for (ev_class = 0; ev_class < evNr; ++ev_class)
803 if (opt2parg_bSet(evOptions[ev_class], NPA, pa))
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)
809 /*if apropriate get list of stepsizes for these eigenvectors*/
810 if (opt2parg_bSet(evStepOptions[ev_class], NPA, pa))
812 evStepList[ev_class] =
813 scan_vecparams(opt2parg_str(evStepOptions[ev_class], NPA, pa), evStepOptions[ev_class], nvecs);
815 else /*if list is not given fill with zeros */
817 snew(evStepList[ev_class], nvecs);
818 for (i = 0; i < nvecs; i++)
820 evStepList[ev_class][i] = 0.0;
824 else if (ev_class == evRADFIX && opt2parg_bSet(evStepOptions[ev_class], NPA, pa))
826 snew(evStepList[ev_class], nvecs);
827 for (i = 0; i < nvecs; i++)
829 evStepList[ev_class][i] = radfix;
832 else if (ev_class == evFLOOD)
834 snew(evStepList[ev_class], nvecs);
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))
840 evStepList[ev_class] = scan_vecparams(opt2parg_str("-constF", NPA, pa), "-constF", nvecs);
845 }; /*to avoid ambiguity */
847 else /* if there are no eigenvectors for this option set list to zero */
849 listen[ev_class] = NULL;
850 snew(listen[ev_class], 1);
851 listen[ev_class][0] = 0;
855 /* print the interpreted list of eigenvectors - to give some feedback*/
856 for (ev_class = 0; ev_class < evNr; ++ev_class)
858 printf("Eigenvector list %7s consists of the indices: ", evOptions[ev_class]);
860 while (listen[ev_class][i])
862 printf("%d ", listen[ev_class][i++]);
868 EigvecFile = opt2fn("-f", NFILE, fnm);
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);
874 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm),
875 title, &top, &ePBC, &xtop, NULL, topbox, 0);
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 */
883 gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
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");
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");
903 get_index(atoms, indexfile, 1, &nfit, &ifit, &grpname);
905 for (i = 0; i < nfit; i++)
907 copy_rvec(xtop[ifit[i]], xref1[i]);
916 if (opt2parg_bSet("-constF", NPA, pa))
918 /* Constant force flooding is special: Most of the normal flooding
919 * options are not needed. */
920 edi_params.flood.bConstForce = TRUE;
924 /* For normal flooding read eigenvalues and store them in evSteplist[evFLOOD] */
926 if (listen[evFLOOD][0] != 0)
928 read_eigenvalues(listen[evFLOOD], opt2fn("-eig", NFILE, fnm), evStepList[evFLOOD], bHesse, kB*T);
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;
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);
945 edi_params.flood.constEfl = constEfl;
946 edi_params.flood.alpha2 = sqr(alpha);
950 edi_params.ned = nav;
952 /*number of system atoms */
953 edi_params.nini = atoms->nr;
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);
961 /* Store target positions in edi_params */
962 if (opt2bSet("-tar", NFILE, fnm))
964 if (0 != listen[evFLOOD][0])
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");
969 get_structure(atoms, indexfile, TargetFile, &edi_params.star, nfit, ifit, nav, index);
973 make_t_edx(&edi_params.star, 0, NULL, index);
976 /* Store origin positions */
977 if (opt2bSet("-ori", NFILE, fnm))
979 get_structure(atoms, indexfile, OriginFile, &edi_params.sori, nfit, ifit, nav, index);
983 make_t_edx(&edi_params.sori, 0, NULL, index);
987 write_the_whole_thing(ffopen(EdiFile, "w"), &edi_params, eigvec1, nvec1, listen, evStepList);