Remove gmx custom fixed int (e.g. gmx_int64_t) types
[gromacs.git] / src / gromacs / gmxana / gmx_make_edi.cpp
blob31c7e6e11f4122a0113ad54511a42a7784d7f4ab
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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/math/functions.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 static 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 static 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 static 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 = nullptr; /*holds the string of the number behind a ','*/
152 char *end = nullptr; /*holds the string of the number behind a '-' */
154 int nvecs = 0; /* counts the number of vectors in the list*/
156 step = nullptr;
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 = nullptr;
166 while ((c = *pos) != 0)
168 switch (status)
170 /* expect a number */
171 case sBefore:
172 if (std::isdigit(c))
174 start = pos;
175 status = sNumber;
176 break;
178 else
180 status = sError;
182 break;
184 /* have read a number, expect ',' or '-' */
185 case sNumber:
186 if (c == ',')
188 /*store number*/
189 srenew(*list, nvecs+1);
190 (*list)[nvecs++] = number = std::strtol(start, nullptr, 10);
191 status = sBefore;
192 if (number == 0)
194 status = sZero;
196 break;
198 else if (c == '-')
200 status = sMinus; break;
202 else if (std::isdigit(c))
204 break;
206 else
208 status = sError;
210 break;
212 /* have read a '-' -> expect a number */
213 case sMinus:
214 if (std::isdigit(c))
216 end = pos;
217 status = sRange; break;
219 else
221 status = sError;
223 break;
225 case sSteppedRange:
226 if (std::isdigit(c))
228 if (step)
230 status = sError; break;
232 else
234 step = pos;
236 status = sRange;
237 break;
239 else
241 status = sError;
243 break;
245 /* have read the number after a minus, expect ',' or ':' */
246 case sRange:
247 if (c == ',')
249 /*store numbers*/
250 end_number = std::strtol(end, nullptr, 10);
251 number = std::strtol(start, nullptr, 10);
252 status = sBefore;
253 if (number == 0)
255 status = sZero; break;
257 if (end_number <= number)
259 status = sSmaller; break;
261 srenew(*list, nvecs+end_number-number+1);
262 if (step)
264 istep = strtol(step, nullptr, 10);
265 step = nullptr;
267 else
269 istep = 1;
271 for (i = number; i <= end_number; i += istep)
273 (*list)[nvecs++] = i;
275 break;
277 else if (c == ':')
279 status = sSteppedRange;
280 break;
282 else if (std::isdigit(c))
284 break;
286 else
288 status = sError;
290 break;
292 /* format error occured */
293 case sError:
294 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %ld with char %c", listname, pos-startpos, *(pos-1));
295 /* logical error occured */
296 case sZero:
297 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %ld: eigenvector 0 is not valid", listname, pos-startpos);
298 case sSmaller:
299 gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %ld: second index %d is not bigger than %d", listname, pos-startpos, end_number, number);
301 ++pos; /* read next character */
302 } /*scanner has finished */
304 /* append zero to list of eigenvectors */
305 srenew(*list, nvecs+1);
306 (*list)[nvecs] = 0;
307 sfree(startpos);
308 return nvecs;
309 } /*sscan_list*/
311 static void write_eigvec(FILE* fp, int natoms, int eig_list[], rvec** eigvecs, int nvec, const char *grouptitle, real steps[])
313 /* eig_list is a zero-terminated list of indices into the eigvecs array.
314 eigvecs are coordinates of eigenvectors
315 grouptitle to write in the comment line
316 steps -- array with stepsizes for evLINFIX, evLINACC and evRADACC
319 int n = 0, i; rvec x;
320 while (eig_list[n++])
322 ; /*count selected eigenvecs*/
325 fprintf(fp, "# NUMBER OF EIGENVECTORS + %s\n %d\n", grouptitle, n-1);
327 /* write list of eigenvector indicess */
328 for (n = 0; eig_list[n]; n++)
330 if (steps)
332 fprintf(fp, "%8d %g\n", eig_list[n], steps[n]);
334 else
336 fprintf(fp, "%8d %g\n", eig_list[n], 1.0);
339 n = 0;
341 /* dump coordinates of the selected eigenvectors */
342 while (eig_list[n])
344 for (i = 0; i < natoms; i++)
346 if (eig_list[n] > nvec)
348 gmx_fatal(FARGS, "Selected eigenvector %d is higher than maximum number %d of available eigenvectors", eig_list[n], nvec);
350 copy_rvec(eigvecs[eig_list[n]-1][i], x);
351 fprintf(fp, "%8.5f %8.5f %8.5f\n", x[XX], x[YY], x[ZZ]);
353 n++;
358 /*enum referring to the different lists of eigenvectors*/
359 enum {
360 evLINFIX, evLINACC, evFLOOD, evRADFIX, evRADACC, evRADCON, evMON, evNr
362 #define oldMAGIC 666
363 #define MAGIC 670
366 static void write_the_whole_thing(FILE* fp, t_edipar *edpars, rvec** eigvecs,
367 int nvec, int *eig_listen[], real* evStepList[])
369 /* write edi-file */
371 /*Header*/
372 fprintf(fp, "#MAGIC\n %d \n#NINI\n %d\n#FITMAS\n %d\n#ANALYSIS_MAS\n %d\n",
373 MAGIC, edpars->nini, int(edpars->fitmas), int(edpars->pcamas));
374 fprintf(fp, "#OUTFRQ\n %d\n#MAXLEN\n %d\n#SLOPECRIT\n %f\n",
375 edpars->outfrq, edpars->maxedsteps, edpars->slope);
376 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",
377 edpars->presteps, edpars->flood.deltaF0, edpars->flood.deltaF, edpars->flood.tau, edpars->flood.constEfl,
378 edpars->flood.alpha2, edpars->flood.kT, int(edpars->flood.bHarmonic), int(edpars->flood.bConstForce));
380 /* Average and reference positions */
381 write_t_edx(fp, edpars->sref, "NREF, XREF");
382 write_t_edx(fp, edpars->sav, "NAV, XAV");
384 /*Eigenvectors */
386 write_eigvec(fp, edpars->ned, eig_listen[evMON], eigvecs, nvec, "COMPONENTS GROUP 1", nullptr);
387 write_eigvec(fp, edpars->ned, eig_listen[evLINFIX], eigvecs, nvec, "COMPONENTS GROUP 2", evStepList[evLINFIX]);
388 write_eigvec(fp, edpars->ned, eig_listen[evLINACC], eigvecs, nvec, "COMPONENTS GROUP 3", evStepList[evLINACC]);
389 write_eigvec(fp, edpars->ned, eig_listen[evRADFIX], eigvecs, nvec, "COMPONENTS GROUP 4", evStepList[evRADFIX]);
390 write_eigvec(fp, edpars->ned, eig_listen[evRADACC], eigvecs, nvec, "COMPONENTS GROUP 5", nullptr);
391 write_eigvec(fp, edpars->ned, eig_listen[evRADCON], eigvecs, nvec, "COMPONENTS GROUP 6", nullptr);
392 write_eigvec(fp, edpars->ned, eig_listen[evFLOOD], eigvecs, nvec, "COMPONENTS GROUP 7", evStepList[evFLOOD]);
395 /*Target and Origin positions */
396 write_t_edx(fp, edpars->star, "NTARGET, XTARGET");
397 write_t_edx(fp, edpars->sori, "NORIGIN, XORIGIN");
400 static int read_conffile(const char *confin, rvec **x)
402 t_topology top;
403 matrix box;
404 printf("read coordnumber from file %s\n", confin);
405 read_tps_conf(confin, &top, nullptr, x, nullptr, box, FALSE);
406 printf("number of coordinates in file %d\n", top.atoms.nr);
407 return top.atoms.nr;
411 static void read_eigenvalues(const int vecs[], const char *eigfile, real values[],
412 gmx_bool bHesse, real kT, int natoms_average_struct)
414 int neig, nrow, i;
415 double **eigval;
417 neig = read_xvg(eigfile, &eigval, &nrow);
419 fprintf(stderr, "Read %d eigenvalues\n", neig);
420 for (i = bHesse ? 6 : 0; i < neig; i++)
422 if (eigval[1][i] < -0.001 && bHesse)
424 fprintf(stderr,
425 "WARNING: The Hessian Matrix has negative eigenvalue %f, we set it to zero (no flooding in this direction)\n\n", eigval[1][i]);
428 if (eigval[1][i] < 0)
430 eigval[1][i] = 0;
433 if (bHesse)
435 for (i = 0; vecs[i]; i++)
437 if (vecs[i] < 7)
439 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");
441 values[i] = eigval[1][vecs[i]-1]/kT;
444 else
446 for (i = 0; vecs[i]; i++)
448 /* Make sure this eigenvalue does not correspond to one of the last 6 eigenvectors of the
449 * covariance matrix. These correspond to the rotational and translational degrees of
450 * freedom and will be zero within numerical accuracy.
452 * Note that the total number of eigenvectors produced by gmx covar depends on:
453 * 1) the total number of degrees of freedom of the system (3N, with N the number of atoms)
454 * 2) the number S of independent configurations fed into gmx covar.
455 * For long trajectories with lots of frames, usually S >= 3N + 1, so that one indeed gets
456 * 3N eigenvalues (of which the last 6 will have zero eigenvalues).
457 * For S < 3N + 1, however, the covariance matrix becomes rank deficient, and the number
458 * of possible eigenvalues is just S - 1. Since in make_edi we only know N but not S, we can
459 * only warn the user if he picked one of the last 6 of 3N.
461 if (vecs[i] > ( 3 * natoms_average_struct - 6 ))
463 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");
465 values[i] = 1/eigval[1][vecs[i]-1];
468 /* free memory */
469 for (i = 0; i < nrow; i++)
471 sfree(eigval[i]);
473 sfree(eigval);
477 static real *scan_vecparams(const char *str, const char * par, int nvecs)
479 char f0[256], f1[256]; /*format strings adapted every pass of the loop*/
480 double d;
481 int i;
482 real *vec_params;
484 snew(vec_params, nvecs);
485 if (str)
487 f0[0] = '\0';
488 for (i = 0; (i < nvecs); i++)
490 std::strcpy(f1, f0); /*f0 is the format string for the "to-be-ignored" numbers*/
491 std::strcat(f1, "%lf"); /*and f1 to read the actual number in this pass of the loop*/
492 if (sscanf(str, f1, &d) != 1)
494 gmx_fatal(FARGS, "Not enough elements for %s parameter (I need %d)", par, nvecs);
496 vec_params[i] = d;
497 std::strcat(f0, "%*s");
500 return vec_params;
504 static void init_edx(struct edix *edx)
506 edx->nr = 0;
507 snew(edx->x, 1);
508 snew(edx->anrs, 1);
511 static void filter2edx(struct edix *edx, int nindex, int index[], int ngro,
512 const int igro[], const rvec *x, const char* structure)
514 /* filter2edx copies coordinates from x to edx which are given in index
517 int pos, i;
518 int ix = edx->nr;
519 edx->nr += nindex;
520 srenew(edx->x, edx->nr);
521 srenew(edx->anrs, edx->nr);
522 for (i = 0; i < nindex; i++, ix++)
524 for (pos = 0; pos < ngro-1 && igro[pos] != index[i]; ++pos)
527 ; /*search element in igro*/
528 if (igro[pos] != index[i])
530 gmx_fatal(FARGS, "Couldn't find atom with index %d in structure %s", index[i], structure);
532 edx->anrs[ix] = index[i];
533 copy_rvec(x[pos], edx->x[ix]);
537 static void get_structure(const t_atoms *atoms, const char *IndexFile,
538 const char *StructureFile, struct edix *edx, int nfit,
539 int ifit[], int nav, int index[])
541 int *igro; /*index corresponding to target or origin structure*/
542 int ngro;
543 int ntar;
544 rvec *xtar;
545 char * grpname;
548 ntar = read_conffile(StructureFile, &xtar);
549 printf("Select an index group of %d elements that corresponds to the atoms in the structure file %s\n",
550 ntar, StructureFile);
551 get_index(atoms, IndexFile, 1, &ngro, &igro, &grpname);
552 if (ngro != ntar)
554 gmx_fatal(FARGS, "You selected an index group with %d elements instead of %d", ngro, ntar);
556 init_edx(edx);
557 filter2edx(edx, nfit, ifit, ngro, igro, xtar, StructureFile);
559 /* If average and reference/fitting structure differ, append the average structure as well */
560 if (ifit != index) /*if fit structure is different append these coordinates, too -- don't mind duplicates*/
562 filter2edx(edx, nav, index, ngro, igro, xtar, StructureFile);
566 int gmx_make_edi(int argc, char *argv[])
569 static const char *desc[] = {
570 "[THISMODULE] generates an essential dynamics (ED) sampling input file to be used with [TT]mdrun[tt]",
571 "based on eigenvectors of a covariance matrix ([gmx-covar]) or from a",
572 "normal modes analysis ([gmx-nmeig]).",
573 "ED sampling can be used to manipulate the position along collective coordinates",
574 "(eigenvectors) of (biological) macromolecules during a simulation. Particularly,",
575 "it may be used to enhance the sampling efficiency of MD simulations by stimulating",
576 "the system to explore new regions along these collective coordinates. A number",
577 "of different algorithms are implemented to drive the system along the eigenvectors",
578 "([TT]-linfix[tt], [TT]-linacc[tt], [TT]-radfix[tt], [TT]-radacc[tt], [TT]-radcon[tt]),",
579 "to keep the position along a certain (set of) coordinate(s) fixed ([TT]-linfix[tt]),",
580 "or to only monitor the projections of the positions onto",
581 "these coordinates ([TT]-mon[tt]).[PAR]",
582 "References:[PAR]",
583 "A. Amadei, A.B.M. Linssen, B.L. de Groot, D.M.F. van Aalten and ",
584 "H.J.C. Berendsen; An efficient method for sampling the essential subspace ",
585 "of proteins., J. Biomol. Struct. Dyn. 13:615-626 (1996)[PAR]",
586 "B.L. de Groot, A. Amadei, D.M.F. van Aalten and H.J.C. Berendsen; ",
587 "Towards an exhaustive sampling of the configurational spaces of the ",
588 "two forms of the peptide hormone guanylin,",
589 "J. Biomol. Struct. Dyn. 13 : 741-751 (1996)[PAR]",
590 "B.L. de Groot, A.Amadei, R.M. Scheek, N.A.J. van Nuland and H.J.C. Berendsen; ",
591 "An extended sampling of the configurational space of HPr from E. coli",
592 "Proteins: Struct. Funct. Gen. 26: 314-322 (1996)",
593 "[PAR]You will be prompted for one or more index groups that correspond to the eigenvectors,",
594 "reference structure, target positions, etc.[PAR]",
596 "[TT]-mon[tt]: monitor projections of the coordinates onto selected eigenvectors.[PAR]",
597 "[TT]-linfix[tt]: perform fixed-step linear expansion along selected eigenvectors.[PAR]",
598 "[TT]-linacc[tt]: perform acceptance linear expansion along selected eigenvectors.",
599 "(steps in the desired directions will be accepted, others will be rejected).[PAR]",
600 "[TT]-radfix[tt]: perform fixed-step radius expansion along selected eigenvectors.[PAR]",
601 "[TT]-radacc[tt]: perform acceptance radius expansion along selected eigenvectors.",
602 "(steps in the desired direction will be accepted, others will be rejected).",
603 "[BB]Note:[bb] by default the starting MD structure will be taken as origin of the first",
604 "expansion cycle for radius expansion. If [TT]-ori[tt] is specified, you will be able",
605 "to read in a structure file that defines an external origin.[PAR]",
606 "[TT]-radcon[tt]: perform acceptance radius contraction along selected eigenvectors",
607 "towards a target structure specified with [TT]-tar[tt].[PAR]",
608 "NOTE: each eigenvector can be selected only once. [PAR]",
609 "[TT]-outfrq[tt]: frequency (in steps) of writing out projections etc. to [REF].xvg[ref] file[PAR]",
610 "[TT]-slope[tt]: minimal slope in acceptance radius expansion. A new expansion",
611 "cycle will be started if the spontaneous increase of the radius (in nm/step)",
612 "is less than the value specified.[PAR]",
613 "[TT]-maxedsteps[tt]: maximum number of steps per cycle in radius expansion",
614 "before a new cycle is started.[PAR]",
615 "Note on the parallel implementation: since ED sampling is a 'global' thing",
616 "(collective coordinates etc.), at least on the 'protein' side, ED sampling",
617 "is not very parallel-friendly from an implementation point of view. Because",
618 "parallel ED requires some extra communication, expect the performance to be",
619 "lower as in a free MD simulation, especially on a large number of ranks and/or",
620 "when the ED group contains a lot of atoms. [PAR]",
621 "Please also note that if your ED group contains more than a single protein,",
622 "then the [REF].tpr[ref] file must contain the correct PBC representation of the ED group.",
623 "Take a look on the initial RMSD from the reference structure, which is printed",
624 "out at the start of the simulation; if this is much higher than expected, one",
625 "of the ED molecules might be shifted by a box vector. [PAR]",
626 "All ED-related output of [TT]mdrun[tt] (specify with [TT]-eo[tt]) is written to a [REF].xvg[ref] file",
627 "as a function of time in intervals of OUTFRQ steps.[PAR]",
628 "[BB]Note[bb] that you can impose multiple ED constraints and flooding potentials in",
629 "a single simulation (on different molecules) if several [REF].edi[ref] files were concatenated",
630 "first. The constraints are applied in the order they appear in the [REF].edi[ref] file. ",
631 "Depending on what was specified in the [REF].edi[ref] input file, the output file contains for each ED dataset",
633 " * the RMSD of the fitted molecule to the reference structure (for atoms involved in fitting prior to calculating the ED constraints)",
634 " * projections of the positions onto selected eigenvectors",
636 "FLOODING:[PAR]",
637 "with [TT]-flood[tt], you can specify which eigenvectors are used to compute a flooding potential,",
638 "which will lead to extra forces expelling the structure out of the region described",
639 "by the covariance matrix. If you switch -restrain the potential is inverted and the structure",
640 "is kept in that region.",
641 "[PAR]",
642 "The origin is normally the average structure stored in the [TT]eigvec.trr[tt] file.",
643 "It can be changed with [TT]-ori[tt] to an arbitrary position in configuration space.",
644 "With [TT]-tau[tt], [TT]-deltaF0[tt], and [TT]-Eflnull[tt] you control the flooding behaviour.",
645 "Efl is the flooding strength, it is updated according to the rule of adaptive flooding.",
646 "Tau is the time constant of adaptive flooding, high [GRK]tau[grk] means slow adaption (i.e. growth). ",
647 "DeltaF0 is the flooding strength you want to reach after tau ps of simulation.",
648 "To use constant Efl set [TT]-tau[tt] to zero.",
649 "[PAR]",
650 "[TT]-alpha[tt] is a fudge parameter to control the width of the flooding potential. A value of 2 has been found",
651 "to give good results for most standard cases in flooding of proteins.",
652 "[GRK]alpha[grk] basically accounts for incomplete sampling, if you sampled further the width of the ensemble would",
653 "increase, this is mimicked by [GRK]alpha[grk] > 1.",
654 "For restraining, [GRK]alpha[grk] < 1 can give you smaller width in the restraining potential.",
655 "[PAR]",
656 "RESTART and FLOODING:",
657 "If you want to restart a crashed flooding simulation please find the values deltaF and Efl in",
658 "the output file and manually put them into the [REF].edi[ref] file under DELTA_F0 and EFL_NULL."
661 /* Save all the params in this struct and then save it in an edi file.
662 * ignoring fields nmass,massnrs,mass,tmass,nfit,fitnrs,edo
664 static t_edipar edi_params;
666 enum {
667 evStepNr = evRADFIX + 1
669 static const char* evSelections[evNr] = {nullptr, nullptr, nullptr, nullptr, nullptr, nullptr};
670 static const char* evOptions[evNr] = {"-linfix", "-linacc", "-flood", "-radfix", "-radacc", "-radcon", "-mon"};
671 static const char* evParams[evStepNr] = {nullptr, nullptr};
672 static const char* evStepOptions[evStepNr] = {"-linstep", "-accdir", "-not_used", "-radstep"};
673 static const char* ConstForceStr;
674 static real * evStepList[evStepNr];
675 static real radstep = 0.0;
676 static real deltaF0 = 150;
677 static real deltaF = 0;
678 static real tau = .1;
679 static real constEfl = 0.0;
680 static real alpha = 1;
681 static int eqSteps = 0;
682 static int * listen[evNr];
683 static real T = 300.0;
684 const real kB = 2.5 / 300.0; /* k_boltzmann in MD units */
685 static gmx_bool bRestrain = FALSE;
686 static gmx_bool bHesse = FALSE;
687 static gmx_bool bHarmonic = FALSE;
688 t_pargs pa[] = {
689 { "-mon", FALSE, etSTR, {&evSelections[evMON]},
690 "Indices of eigenvectors for projections of x (e.g. 1,2-5,9) or 1-100:10 means 1 11 21 31 ... 91" },
691 { "-linfix", FALSE, etSTR, {&evSelections[0]},
692 "Indices of eigenvectors for fixed increment linear sampling" },
693 { "-linacc", FALSE, etSTR, {&evSelections[1]},
694 "Indices of eigenvectors for acceptance linear sampling" },
695 { "-radfix", FALSE, etSTR, {&evSelections[3]},
696 "Indices of eigenvectors for fixed increment radius expansion" },
697 { "-radacc", FALSE, etSTR, {&evSelections[4]},
698 "Indices of eigenvectors for acceptance radius expansion" },
699 { "-radcon", FALSE, etSTR, {&evSelections[5]},
700 "Indices of eigenvectors for acceptance radius contraction" },
701 { "-flood", FALSE, etSTR, {&evSelections[2]},
702 "Indices of eigenvectors for flooding"},
703 { "-outfrq", FALSE, etINT, {&edi_params.outfrq},
704 "Frequency (in steps) of writing output in [REF].xvg[ref] file" },
705 { "-slope", FALSE, etREAL, { &edi_params.slope},
706 "Minimal slope in acceptance radius expansion"},
707 { "-linstep", FALSE, etSTR, {&evParams[0]},
708 "Stepsizes (nm/step) for fixed increment linear sampling (put in quotes! \"1.0 2.3 5.1 -3.1\")"},
709 { "-accdir", FALSE, etSTR, {&evParams[1]},
710 "Directions for acceptance linear sampling - only sign counts! (put in quotes! \"-1 +1 -1.1\")"},
711 { "-radstep", FALSE, etREAL, {&radstep},
712 "Stepsize (nm/step) for fixed increment radius expansion"},
713 { "-maxedsteps", FALSE, etINT, {&edi_params.maxedsteps},
714 "Maximum number of steps per cycle" },
715 { "-eqsteps", FALSE, etINT, {&eqSteps},
716 "Number of steps to run without any perturbations "},
717 { "-deltaF0", FALSE, etREAL, {&deltaF0},
718 "Target destabilization energy for flooding"},
719 { "-deltaF", FALSE, etREAL, {&deltaF},
720 "Start deltaF with this parameter - default 0, nonzero values only needed for restart"},
721 { "-tau", FALSE, etREAL, {&tau},
722 "Coupling constant for adaption of flooding strength according to deltaF0, 0 = infinity i.e. constant flooding strength"},
723 { "-Eflnull", FALSE, etREAL, {&constEfl},
724 "The starting value of the flooding strength. The flooding strength is updated "
725 "according to the adaptive flooding scheme. For a constant flooding strength use [TT]-tau[tt] 0. "},
726 { "-T", FALSE, etREAL, {&T},
727 "T is temperature, the value is needed if you want to do flooding "},
728 { "-alpha", FALSE, etREAL, {&alpha},
729 "Scale width of gaussian flooding potential with alpha^2 "},
730 { "-restrain", FALSE, etBOOL, {&bRestrain},
731 "Use the flooding potential with inverted sign -> effects as quasiharmonic restraining potential"},
732 { "-hessian", FALSE, etBOOL, {&bHesse},
733 "The eigenvectors and eigenvalues are from a Hessian matrix"},
734 { "-harmonic", FALSE, etBOOL, {&bHarmonic},
735 "The eigenvalues are interpreted as spring constant"},
736 { "-constF", FALSE, etSTR, {&ConstForceStr},
737 "Constant force flooding: manually set the forces for the eigenvectors selected with -flood "
738 "(put in quotes! \"1.0 2.3 5.1 -3.1\"). No other flooding parameters are needed when specifying the forces directly."}
740 #define NPA asize(pa)
742 rvec *xref1;
743 int nvec1, *eignr1 = nullptr;
744 rvec *xav1, **eigvec1 = nullptr;
745 t_atoms *atoms = nullptr;
746 int nav; /* Number of atoms in the average structure */
747 char *grpname;
748 const char *indexfile;
749 int i;
750 int *index, *ifit;
751 int nfit; /* Number of atoms in the reference/fit structure */
752 int ev_class; /* parameter _class i.e. evMON, evRADFIX etc. */
753 int nvecs;
754 real *eigval1 = nullptr; /* in V3.3 this is parameter of read_eigenvectors */
756 const char *EdiFile;
757 const char *TargetFile;
758 const char *OriginFile;
759 const char *EigvecFile;
761 gmx_output_env_t *oenv;
763 /*to read topology file*/
764 t_topology top;
765 int ePBC;
766 matrix topbox;
767 rvec *xtop;
768 gmx_bool bFit1;
770 t_filenm fnm[] = {
771 { efTRN, "-f", "eigenvec", ffREAD },
772 { efXVG, "-eig", "eigenval", ffOPTRD },
773 { efTPS, nullptr, nullptr, ffREAD },
774 { efNDX, nullptr, nullptr, ffOPTRD },
775 { efSTX, "-tar", "target", ffOPTRD},
776 { efSTX, "-ori", "origin", ffOPTRD},
777 { efEDI, "-o", "sam", ffWRITE }
779 #define NFILE asize(fnm)
780 edi_params.outfrq = 100; edi_params.slope = 0.0; edi_params.maxedsteps = 0;
781 if (!parse_common_args(&argc, argv, 0,
782 NFILE, fnm, NPA, pa, asize(desc), desc, 0, nullptr, &oenv))
784 return 0;
787 indexfile = ftp2fn_null(efNDX, NFILE, fnm);
788 EdiFile = ftp2fn(efEDI, NFILE, fnm);
789 TargetFile = opt2fn_null("-tar", NFILE, fnm);
790 OriginFile = opt2fn_null("-ori", NFILE, fnm);
793 for (ev_class = 0; ev_class < evNr; ++ev_class)
795 if (opt2parg_bSet(evOptions[ev_class], NPA, pa))
797 /*get list of eigenvectors*/
798 nvecs = sscan_list(&(listen[ev_class]), opt2parg_str(evOptions[ev_class], NPA, pa), evOptions[ev_class]);
799 if (ev_class < evStepNr-2)
801 /*if apropriate get list of stepsizes for these eigenvectors*/
802 if (opt2parg_bSet(evStepOptions[ev_class], NPA, pa))
804 evStepList[ev_class] =
805 scan_vecparams(opt2parg_str(evStepOptions[ev_class], NPA, pa), evStepOptions[ev_class], nvecs);
807 else /*if list is not given fill with zeros */
809 snew(evStepList[ev_class], nvecs);
810 for (i = 0; i < nvecs; i++)
812 evStepList[ev_class][i] = 0.0;
816 else if (ev_class == evRADFIX)
818 snew(evStepList[ev_class], nvecs);
819 for (i = 0; i < nvecs; i++)
821 evStepList[ev_class][i] = radstep;
824 else if (ev_class == evFLOOD)
826 snew(evStepList[ev_class], nvecs);
828 /* Are we doing constant force flooding? In that case, we read in
829 * the fproj values from the command line */
830 if (opt2parg_bSet("-constF", NPA, pa))
832 evStepList[ev_class] = scan_vecparams(opt2parg_str("-constF", NPA, pa), "-constF", nvecs);
835 else
837 }; /*to avoid ambiguity */
839 else /* if there are no eigenvectors for this option set list to zero */
841 listen[ev_class] = nullptr;
842 snew(listen[ev_class], 1);
843 listen[ev_class][0] = 0;
847 /* print the interpreted list of eigenvectors - to give some feedback*/
848 for (ev_class = 0; ev_class < evNr; ++ev_class)
850 printf("Eigenvector list %7s consists of the indices: ", evOptions[ev_class]);
851 i = 0;
852 while (listen[ev_class][i])
854 printf("%d ", listen[ev_class][i++]);
856 printf("\n");
859 EigvecFile = opt2fn("-f", NFILE, fnm);
861 /*read eigenvectors from eigvec.trr*/
862 read_eigenvectors(EigvecFile, &nav, &bFit1,
863 &xref1, &edi_params.fitmas, &xav1, &edi_params.pcamas, &nvec1, &eignr1, &eigvec1, &eigval1);
865 read_tps_conf(ftp2fn(efTPS, NFILE, fnm),
866 &top, &ePBC, &xtop, nullptr, topbox, 0);
867 atoms = &top.atoms;
870 printf("\nSelect an index group of %d elements that corresponds to the eigenvectors\n", nav);
871 get_index(atoms, indexfile, 1, &i, &index, &grpname); /*if indexfile != NULL parameter 'atoms' is ignored */
872 if (i != nav)
874 gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
875 i, nav);
877 printf("\n");
880 if (xref1 == nullptr)
882 if (bFit1)
884 /* if g_covar used different coordinate groups to fit and to do the PCA */
885 printf("\nNote: the structure in %s should be the same\n"
886 " as the one used for the fit in g_covar\n", ftp2fn(efTPS, NFILE, fnm));
887 printf("\nSelect the index group that was used for the least squares fit in g_covar\n");
889 else
891 printf("\nNote: Apparently no fitting was done in g_covar.\n"
892 " However, you need to select a reference group for fitting in mdrun\n");
894 get_index(atoms, indexfile, 1, &nfit, &ifit, &grpname);
895 snew(xref1, nfit);
896 for (i = 0; i < nfit; i++)
898 copy_rvec(xtop[ifit[i]], xref1[i]);
901 else
903 nfit = nav;
904 ifit = index;
907 if (opt2parg_bSet("-constF", NPA, pa))
909 /* Constant force flooding is special: Most of the normal flooding
910 * options are not needed. */
911 edi_params.flood.bConstForce = TRUE;
913 else
915 /* For normal flooding read eigenvalues and store them in evSteplist[evFLOOD] */
917 if (listen[evFLOOD][0] != 0)
919 read_eigenvalues(listen[evFLOOD], opt2fn("-eig", NFILE, fnm), evStepList[evFLOOD], bHesse, kB*T, nav);
922 edi_params.flood.tau = tau;
923 edi_params.flood.deltaF0 = deltaF0;
924 edi_params.flood.deltaF = deltaF;
925 edi_params.presteps = eqSteps;
926 edi_params.flood.kT = kB*T;
927 edi_params.flood.bHarmonic = bHarmonic;
928 if (bRestrain)
930 /* Trick: invert sign of Efl and alpha2 then this will give the same sign in the exponential and inverted sign outside */
931 edi_params.flood.constEfl = -constEfl;
932 edi_params.flood.alpha2 = -gmx::square(alpha);
934 else
936 edi_params.flood.constEfl = constEfl;
937 edi_params.flood.alpha2 = gmx::square(alpha);
941 edi_params.ned = nav;
943 /*number of system atoms */
944 edi_params.nini = atoms->nr;
947 /*store reference and average structure in edi_params*/
948 make_t_edx(&edi_params.sref, nfit, xref1, ifit );
949 make_t_edx(&edi_params.sav, nav, xav1, index);
952 /* Store target positions in edi_params */
953 if (opt2bSet("-tar", NFILE, fnm))
955 if (0 != listen[evFLOOD][0])
957 fprintf(stderr, "\nNote: Providing a TARGET structure has no effect when using flooding.\n"
958 " You may want to use -ori to define the flooding potential center.\n\n");
960 get_structure(atoms, indexfile, TargetFile, &edi_params.star, nfit, ifit, nav, index);
962 else
964 make_t_edx(&edi_params.star, 0, nullptr, index);
967 /* Store origin positions */
968 if (opt2bSet("-ori", NFILE, fnm))
970 get_structure(atoms, indexfile, OriginFile, &edi_params.sori, nfit, ifit, nav, index);
972 else
974 make_t_edx(&edi_params.sori, 0, nullptr, index);
977 /* Write edi-file */
978 write_the_whole_thing(gmx_ffopen(EdiFile, "w"), &edi_params, eigvec1, nvec1, listen, evStepList);
980 return 0;