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