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