3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
11 * The make_edi program was generously contributed by Oliver Lange, based
12 * on the code from g_anaeig. You can reach him as olange@gwdg.de.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
32 * Gromacs Runs One Microsecond At Cannonball Speeds
48 #include "gmx_fatal.h"
67 /* Suppress Cygwin compiler warnings from using newlib version of
77 gmx_bool bConstForce; /* Do constant force flooding instead of
78 evaluating a flooding potential */
87 /* This type is for the average, reference, target, and origin structure */
90 int nr; /* number of atoms this structure contains */
91 int *anrs; /* atom index numbers */
92 rvec *x; /* positions */
93 real *sqrtm; /* sqrt of the masses used for mass-
94 * weighting of analysis */
100 int nini; /* total Nr of atoms */
101 gmx_bool fitmas; /* true if trans fit with cm */
102 gmx_bool pcamas; /* true if mass-weighted PCA */
103 int presteps; /* number of steps to run without any
104 * perturbations ... just monitoring */
105 int outfrq; /* freq (in steps) of writing to edo */
106 int maxedsteps; /* max nr of steps per cycle */
107 struct edix sref; /* reference positions, to these fitting
109 struct edix sav; /* average positions */
110 struct edix star; /* target positions */
111 struct edix sori; /* origin positions */
112 real slope; /* minimal slope in acceptance radexp */
113 int ned; /* Nr of atoms in essdyn buffer */
114 t_edflood flood; /* parameters especially for flooding */
119 void make_t_edx(struct edix *edx, int natoms, rvec *pos, atom_id index[]) {
125 void write_t_edx(FILE *fp, struct edix edx, const char *comment) {
126 /*here we copy only the pointers into the t_edx struct
127 no data is copied and edx.box is ignored */
129 fprintf(fp,"#%s \n %d \n",comment,edx.nr);
130 for (i=0;i<edx.nr;i++) {
131 fprintf(fp,"%d %f %f %f\n",(edx.anrs)[i]+1,(edx.x)[i][XX],(edx.x)[i][YY],(edx.x)[i][ZZ]);
135 int sscan_list(int *list[], const char *str, const char *listname) {
136 /*this routine scans a string of the form 1,3-6,9 and returns the
137 selected numbers (in this case 1 3 4 5 6 9) in NULL-terminated array of integers.
138 memory for this list will be allocated in this routine -- sscan_list expects *list to
141 listname is a string used in the errormessage*/
146 char *pos,*startpos,*step;
149 /*enums to define the different lexical stati */
150 enum { sBefore, sNumber, sMinus, sRange, sZero, sSmaller, sError, sSteppedRange };
152 int status=sBefore; /*status of the deterministic automat to scan str */
156 char *start=NULL; /*holds the string of the number behind a ','*/
157 char *end=NULL; /*holds the string of the number behind a '-' */
159 int nvecs=0; /* counts the number of vectors in the list*/
171 while ((c=*pos)!=0) {
173 /* expect a number */
174 case sBefore: if (isdigit(c)) {
178 } else status=sError; break;
180 /* have read a number, expect ',' or '-' */
181 case sNumber: if (c==',') {
183 srenew(*list,nvecs+1);
184 (*list)[nvecs++]=number=strtol(start,NULL,10);
189 } else if (c=='-') { status=sMinus; break; }
190 else if (isdigit(c)) break;
191 else status=sError; break;
193 /* have read a '-' -> expect a number */
197 status=sRange; break;
198 } else status=sError; break;
203 status=sError; break;
208 } else status=sError; break;
210 /* have read the number after a minus, expect ',' or ':' */
214 end_number=strtol(end,NULL,10);
215 number=strtol(start,NULL,10);
220 if (end_number<=number) {
221 status=sSmaller; break;
223 srenew(*list,nvecs+end_number-number+1);
225 istep=strtol(step,NULL,10);
228 for (i=number;i<=end_number;i+=istep)
232 status = sSteppedRange;
234 } else if (isdigit(c)) break; else status=sError; break;
236 /* format error occured */
238 gmx_fatal(FARGS,"Error in the list of eigenvectors for %s at pos %d with char %c",listname,pos-startpos,*(pos-1));
240 /* logical error occured */
242 gmx_fatal(FARGS,"Error in the list of eigenvectors for %s at pos %d: eigenvector 0 is not valid",listname,pos-startpos);
244 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);
247 ++pos; /* read next character */
248 } /*scanner has finished */
250 /* append zero to list of eigenvectors */
251 srenew(*list,nvecs+1);
257 void write_eigvec(FILE* fp, int natoms, int eig_list[], rvec** eigvecs,int nvec, const char *grouptitle,real steps[]) {
258 /* eig_list is a zero-terminated list of indices into the eigvecs array.
259 eigvecs are coordinates of eigenvectors
260 grouptitle to write in the comment line
261 steps -- array with stepsizes for evLINFIX, evLINACC and evRADACC
266 while (eig_list[n++]); /*count selected eigenvecs*/
268 fprintf(fp,"# NUMBER OF EIGENVECTORS + %s\n %d\n",grouptitle,n-1);
270 /* write list of eigenvector indicess */
271 for(n=0;eig_list[n];n++) {
273 fprintf(fp,"%8d %g\n",eig_list[n],steps[n]);
275 fprintf(fp,"%8d %g\n",eig_list[n],1.0);
279 /* dump coordinates of the selected eigenvectors */
280 while (eig_list[n]) {
282 for (i=0; i<natoms; i++) {
283 if (eig_list[n]>nvec)
284 gmx_fatal(FARGS,"Selected eigenvector %d is higher than maximum number %d of available eigenvectors",eig_list[n],nvec);
285 copy_rvec(eigvecs[eig_list[n]-1][i],x);
287 fprintf(fp,"%8.5f %8.5f %8.5f\n",x[XX],x[YY],x[ZZ]);
294 /*enum referring to the different lists of eigenvectors*/
295 enum { evLINFIX, evLINACC, evFLOOD, evRADFIX, evRADACC, evRADCON , evMON, evNr };
300 void write_the_whole_thing(FILE* fp, t_edipar *edpars, rvec** eigvecs,
301 int nvec, int *eig_listen[], real* evStepList[])
306 fprintf(fp,"#MAGIC\n %d \n#NINI\n %d\n#FITMAS\n %d\n#ANALYSIS_MAS\n %d\n",
307 MAGIC,edpars->nini,edpars->fitmas,edpars->pcamas);
308 fprintf(fp,"#OUTFRQ\n %d\n#MAXLEN\n %d\n#SLOPECRIT\n %f\n",
309 edpars->outfrq,edpars->maxedsteps,edpars->slope);
310 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",
311 edpars->presteps,edpars->flood.deltaF0,edpars->flood.deltaF,edpars->flood.tau,edpars->flood.constEfl,
312 edpars->flood.alpha2,edpars->flood.kT,edpars->flood.bHarmonic,edpars->flood.bConstForce);
314 /* Average and reference positions */
315 write_t_edx(fp,edpars->sref,"NREF, XREF");
316 write_t_edx(fp,edpars->sav,"NAV, XAV");
320 write_eigvec(fp, edpars->ned, eig_listen[evMON],eigvecs,nvec,"COMPONENTS GROUP 1",NULL);
321 write_eigvec(fp, edpars->ned, eig_listen[evLINFIX],eigvecs,nvec,"COMPONENTS GROUP 2",evStepList[evLINFIX]);
322 write_eigvec(fp, edpars->ned, eig_listen[evLINACC],eigvecs,nvec,"COMPONENTS GROUP 3",evStepList[evLINACC]);
323 write_eigvec(fp, edpars->ned, eig_listen[evRADFIX],eigvecs,nvec,"COMPONENTS GROUP 4",evStepList[evRADFIX]);
324 write_eigvec(fp, edpars->ned, eig_listen[evRADACC],eigvecs,nvec,"COMPONENTS GROUP 5",NULL);
325 write_eigvec(fp, edpars->ned, eig_listen[evRADCON],eigvecs,nvec,"COMPONENTS GROUP 6",NULL);
326 write_eigvec(fp, edpars->ned, eig_listen[evFLOOD],eigvecs,nvec,"COMPONENTS GROUP 7",evStepList[evFLOOD]);
329 /*Target and Origin positions */
330 write_t_edx(fp,edpars->star,"NTARGET, XTARGET");
331 write_t_edx(fp,edpars->sori,"NORIGIN, XORIGIN");
334 int read_conffile(const char *confin,char *title,rvec *x[])
336 /* read coordinates out of STX file */
340 printf("read coordnumber from file %s\n",confin);
341 get_stx_coordnum(confin,&natoms);
342 printf("number of coordinates in file %d\n",natoms);
343 /* if (natoms != ncoords)
344 gmx_fatal(FARGS,"number of coordinates in coordinate file (%s, %d)\n"
345 " does not match topology (= %d)",
346 confin,natoms,ncoords);
348 /* make space for coordinates and velocities */
349 init_t_atoms(&confat,natoms,FALSE);
351 read_stx_conf(confin,title,&confat,*x,NULL,NULL,box);
356 void read_eigenvalues(int vecs[],const char *eigfile, real values[],
357 gmx_bool bHesse, real kT)
362 neig = read_xvg(eigfile,&eigval,&nrow);
364 fprintf(stderr,"Read %d eigenvalues\n",neig);
365 for(i=bHesse ? 6 : 0 ; i<neig; i++) {
366 if (eigval[1][i] < -0.001 && bHesse)
368 "WARNING: The Hessian Matrix has negative eigenvalue %f, we set it to zero (no flooding in this direction)\n\n",eigval[1][i]);
370 if (eigval[1][i] < 0)
374 for (i=0; vecs[i]; i++) {
376 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");
377 values[i]=eigval[1][vecs[i]-1]/kT;
380 for (i=0; vecs[i]; i++) {
381 if (vecs[i]>(neig-6))
382 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");
383 values[i]=1/eigval[1][vecs[i]-1];
386 for (i=0; i<nrow; i++)
392 static real *scan_vecparams(const char *str,const char * par, int nvecs)
394 char f0[256],f1[256]; /*format strings adapted every pass of the loop*/
399 snew(vec_params,nvecs);
402 for(i=0; (i<nvecs); i++) {
403 strcpy(f1,f0); /*f0 is the format string for the "to-be-ignored" numbers*/
404 strcat(f1,"%lf"); /*and f1 to read the actual number in this pass of the loop*/
405 if (sscanf(str,f1,&d) != 1)
406 gmx_fatal(FARGS,"Not enough elements for %s parameter (I need %d)",par,nvecs);
416 void init_edx(struct edix *edx) {
422 void filter2edx(struct edix *edx,int nindex, atom_id index[],int ngro,
423 atom_id igro[],rvec *x,const char* structure)
425 /* filter2edx copies coordinates from x to edx which are given in index
431 srenew(edx->x,edx->nr);
432 srenew(edx->anrs,edx->nr);
433 for (i=0;i<nindex;i++,ix++) {
434 for (pos=0; pos<ngro-1 && igro[pos]!=index[i] ; ++pos) {}; /*search element in igro*/
435 if (igro[pos]!=index[i])
436 gmx_fatal(FARGS,"Couldn't find atom with index %d in structure %s",index[i],structure);
437 edx->anrs[ix]=index[i];
438 copy_rvec(x[pos],edx->x[ix]);
442 void get_structure(t_atoms *atoms,const char *IndexFile,
443 const char *StructureFile,struct edix *edx,int nfit,
444 atom_id ifit[],int nav, atom_id index[])
446 atom_id *igro; /*index corresponding to target or origin structure*/
454 ntar=read_conffile(StructureFile,title,&xtar);
455 printf("Select an index group of %d elements that corresponds to the atoms in the structure file %s\n",
457 get_index(atoms,IndexFile,1,&ngro,&igro,&grpname);
459 gmx_fatal(FARGS,"You selected an index group with %d elements instead of %d",ngro,ntar);
461 filter2edx(edx,nfit,ifit,ngro,igro,xtar,StructureFile);
463 /* If average and reference/fitting structure differ, append the average structure as well */
464 if (ifit!=index) /*if fit structure is different append these coordinates, too -- don't mind duplicates*/
465 filter2edx(edx,nav,index,ngro,igro,xtar,StructureFile);
468 int main(int argc,char *argv[])
471 static const char *desc[] = {
472 "[TT]make_edi[tt] generates an essential dynamics (ED) sampling input file to be used with [TT]mdrun[tt]",
473 "based on eigenvectors of a covariance matrix ([TT]g_covar[tt]) or from a",
474 "normal modes anaysis ([TT]g_nmeig[tt]).",
475 "ED sampling can be used to manipulate the position along collective coordinates",
476 "(eigenvectors) of (biological) macromolecules during a simulation. Particularly,",
477 "it may be used to enhance the sampling efficiency of MD simulations by stimulating",
478 "the system to explore new regions along these collective coordinates. A number",
479 "of different algorithms are implemented to drive the system along the eigenvectors",
480 "([TT]-linfix[tt], [TT]-linacc[tt], [TT]-radfix[tt], [TT]-radacc[tt], [TT]-radcon[tt]),",
481 "to keep the position along a certain (set of) coordinate(s) fixed ([TT]-linfix[tt]),",
482 "or to only monitor the projections of the positions onto",
483 "these coordinates ([TT]-mon[tt]).[PAR]",
485 "A. Amadei, A.B.M. Linssen, B.L. de Groot, D.M.F. van Aalten and ",
486 "H.J.C. Berendsen; An efficient method for sampling the essential subspace ",
487 "of proteins., J. Biomol. Struct. Dyn. 13:615-626 (1996)[BR]",
488 "B.L. de Groot, A. Amadei, D.M.F. van Aalten and H.J.C. Berendsen; ",
489 "Towards an exhaustive sampling of the configurational spaces of the ",
490 "two forms of the peptide hormone guanylin,",
491 "J. Biomol. Struct. Dyn. 13 : 741-751 (1996)[BR]",
492 "B.L. de Groot, A.Amadei, R.M. Scheek, N.A.J. van Nuland and H.J.C. Berendsen; ",
493 "An extended sampling of the configurational space of HPr from E. coli",
494 "Proteins: Struct. Funct. Gen. 26: 314-322 (1996)",
495 "[PAR]You will be prompted for one or more index groups that correspond to the eigenvectors,",
496 "reference structure, target positions, etc.[PAR]",
498 "[TT]-mon[tt]: monitor projections of the coordinates onto selected eigenvectors.[PAR]",
499 "[TT]-linfix[tt]: perform fixed-step linear expansion along selected eigenvectors.[PAR]",
500 "[TT]-linacc[tt]: perform acceptance linear expansion along selected eigenvectors.",
501 "(steps in the desired directions will be accepted, others will be rejected).[PAR]",
502 "[TT]-radfix[tt]: perform fixed-step radius expansion along selected eigenvectors.[PAR]",
503 "[TT]-radacc[tt]: perform acceptance radius expansion along selected eigenvectors.",
504 "(steps in the desired direction will be accepted, others will be rejected).",
505 "[BB]Note:[bb] by default the starting MD structure will be taken as origin of the first",
506 "expansion cycle for radius expansion. If [TT]-ori[tt] is specified, you will be able",
507 "to read in a structure file that defines an external origin.[PAR]",
508 "[TT]-radcon[tt]: perform acceptance radius contraction along selected eigenvectors",
509 "towards a target structure specified with [TT]-tar[tt].[PAR]",
510 "NOTE: each eigenvector can be selected only once. [PAR]",
511 "[TT]-outfrq[tt]: frequency (in steps) of writing out projections etc. to [TT].edo[tt] file[PAR]",
512 "[TT]-slope[tt]: minimal slope in acceptance radius expansion. A new expansion",
513 "cycle will be started if the spontaneous increase of the radius (in nm/step)",
514 "is less than the value specified.[PAR]",
515 "[TT]-maxedsteps[tt]: maximum number of steps per cycle in radius expansion",
516 "before a new cycle is started.[PAR]",
517 "Note on the parallel implementation: since ED sampling is a 'global' thing",
518 "(collective coordinates etc.), at least on the 'protein' side, ED sampling",
519 "is not very parallel-friendly from an implentation point of view. Because",
520 "parallel ED requires some extra communication, expect the performance to be",
521 "lower as in a free MD simulation, especially on a large number of nodes. [PAR]",
522 "All output of [TT]mdrun[tt] (specify with [TT]-eo[tt]) is written to a .edo file. In the output",
523 "file, per OUTFRQ step the following information is present: [PAR]",
524 "[TT]*[tt] the step number[BR]",
525 "[TT]*[tt] the number of the ED dataset. ([BB]Note[bb] that you can impose multiple ED constraints in",
526 "a single simulation (on different molecules) if several [TT].edi[tt] files were concatenated",
527 "first. The constraints are applied in the order they appear in the [TT].edi[tt] file.) [BR]",
528 "[TT]*[tt] RMSD (for atoms involved in fitting prior to calculating the ED constraints)[BR]",
529 "* projections of the positions onto selected eigenvectors[BR]",
532 "with [TT]-flood[tt], you can specify which eigenvectors are used to compute a flooding potential,",
533 "which will lead to extra forces expelling the structure out of the region described",
534 "by the covariance matrix. If you switch -restrain the potential is inverted and the structure",
535 "is kept in that region.",
537 "The origin is normally the average structure stored in the [TT]eigvec.trr[tt] file.",
538 "It can be changed with [TT]-ori[tt] to an arbitrary position in configurational space.",
539 "With [TT]-tau[tt], [TT]-deltaF0[tt], and [TT]-Eflnull[tt] you control the flooding behaviour.",
540 "Efl is the flooding strength, it is updated according to the rule of adaptive flooding.",
541 "Tau is the time constant of adaptive flooding, high [GRK]tau[grk] means slow adaption (i.e. growth). ",
542 "DeltaF0 is the flooding strength you want to reach after tau ps of simulation.",
543 "To use constant Efl set [TT]-tau[tt] to zero.",
545 "[TT]-alpha[tt] is a fudge parameter to control the width of the flooding potential. A value of 2 has been found",
546 "to give good results for most standard cases in flooding of proteins.",
547 "[GRK]alpha[grk] basically accounts for incomplete sampling, if you sampled further the width of the ensemble would",
548 "increase, this is mimicked by [GRK]alpha[grk] > 1.",
549 "For restraining, [GRK]alpha[grk] < 1 can give you smaller width in the restraining potential.",
551 "RESTART and FLOODING:",
552 "If you want to restart a crashed flooding simulation please find the values deltaF and Efl in",
553 "the output file and manually put them into the [TT].edi[tt] file under DELTA_F0 and EFL_NULL."
556 /* Save all the params in this struct and then save it in an edi file.
557 * ignoring fields nmass,massnrs,mass,tmass,nfit,fitnrs,edo
559 static t_edipar edi_params;
561 enum { evStepNr = evRADFIX + 1};
562 static const char* evSelections[evNr]= {NULL,NULL,NULL,NULL,NULL,NULL};
563 static const char* evOptions[evNr] = {"-linfix","-linacc","-flood","-radfix","-radacc","-radcon","-mon"};
564 static const char* evParams[evStepNr] ={NULL,NULL};
565 static const char* evStepOptions[evStepNr] = {"-linstep","-accdir","-not_used","-radstep"};
566 static const char* ConstForceStr;
567 static real* evStepList[evStepNr];
568 static real radfix=0.0;
569 static real deltaF0=150;
570 static real deltaF=0;
572 static real constEfl=0.0;
574 static int eqSteps=0;
575 static int* listen[evNr];
577 const real kB = 2.5 / 300.0; /* k_boltzmann in MD units */
578 static gmx_bool bRestrain = FALSE;
579 static gmx_bool bHesse=FALSE;
580 static gmx_bool bHarmonic=FALSE;
582 { "-mon", FALSE, etSTR, {&evSelections[evMON]},
583 "Indices of eigenvectors for projections of x (e.g. 1,2-5,9) or 1-100:10 means 1 11 21 31 ... 91" },
584 { "-linfix", FALSE, etSTR, {&evSelections[0]},
585 "Indices of eigenvectors for fixed increment linear sampling" },
586 { "-linacc", FALSE, etSTR, {&evSelections[1]},
587 "Indices of eigenvectors for acceptance linear sampling" },
588 { "-radfix", FALSE, etSTR, {&evSelections[3]},
589 "Indices of eigenvectors for fixed increment radius expansion" },
590 { "-radacc", FALSE, etSTR, {&evSelections[4]},
591 "Indices of eigenvectors for acceptance radius expansion" },
592 { "-radcon", FALSE, etSTR, {&evSelections[5]},
593 "Indices of eigenvectors for acceptance radius contraction" },
594 { "-flood", FALSE, etSTR, {&evSelections[2]},
595 "Indices of eigenvectors for flooding"},
596 { "-outfrq", FALSE, etINT, {&edi_params.outfrq},
597 "Freqency (in steps) of writing output in [TT].edo[tt] file" },
598 { "-slope", FALSE, etREAL, { &edi_params.slope},
599 "Minimal slope in acceptance radius expansion"},
600 { "-linstep", FALSE, etSTR, {&evParams[0]},
601 "Stepsizes (nm/step) for fixed increment linear sampling (put in quotes! \"1.0 2.3 5.1 -3.1\")"},
602 { "-accdir", FALSE, etSTR, {&evParams[1]},
603 "Directions for acceptance linear sampling - only sign counts! (put in quotes! \"-1 +1 -1.1\")"},
604 { "-radstep", FALSE, etREAL, {&radfix},
605 "Stepsize (nm/step) for fixed increment radius expansion"},
606 { "-maxedsteps", FALSE, etINT, {&edi_params.maxedsteps},
607 "Maximum number of steps per cycle" },
608 { "-eqsteps", FALSE, etINT, {&eqSteps},
609 "Number of steps to run without any perturbations "},
610 { "-deltaF0", FALSE,etREAL, {&deltaF0},
611 "Target destabilization energy for flooding"},
612 { "-deltaF", FALSE, etREAL, {&deltaF},
613 "Start deltaF with this parameter - default 0, nonzero values only needed for restart"},
614 { "-tau", FALSE, etREAL, {&tau},
615 "Coupling constant for adaption of flooding strength according to deltaF0, 0 = infinity i.e. constant flooding strength"},
616 { "-Eflnull", FALSE, etREAL, {&constEfl},
617 "The starting value of the flooding strength. The flooding strength is updated "
618 "according to the adaptive flooding scheme. For a constant flooding strength use [TT]-tau[tt] 0. "},
619 { "-T", FALSE, etREAL, {&T},
620 "T is temperature, the value is needed if you want to do flooding "},
621 { "-alpha",FALSE,etREAL,{&alpha},
622 "Scale width of gaussian flooding potential with alpha^2 "},
623 { "-restrain",FALSE, etBOOL, {&bRestrain},
624 "Use the flooding potential with inverted sign -> effects as quasiharmonic restraining potential"},
625 { "-hessian",FALSE, etBOOL, {&bHesse},
626 "The eigenvectors and eigenvalues are from a Hessian matrix"},
627 { "-harmonic",FALSE, etBOOL, {&bHarmonic},
628 "The eigenvalues are interpreted as spring constant"},
629 { "-constF", FALSE, etSTR, {&ConstForceStr},
630 "Constant force flooding: manually set the forces for the eigenvectors selected with -flood "
631 "(put in quotes! \"1.0 2.3 5.1 -3.1\"). No other flooding parameters are needed when specifying the forces directly."}
633 #define NPA asize(pa)
636 int nvec1,*eignr1=NULL;
637 rvec *xav1,**eigvec1=NULL;
639 int nav; /* Number of atoms in the average structure */
641 const char *indexfile;
643 atom_id *index,*ifit;
644 int nfit; /* Number of atoms in the reference/fit structure */
645 int ev_class; /* parameter _class i.e. evMON, evRADFIX etc. */
647 real *eigval1=NULL; /* in V3.3 this is parameter of read_eigenvectors */
650 const char *TargetFile;
651 const char *OriginFile;
652 const char *EigvecFile;
656 /*to read topology file*/
662 gmx_bool bTop, bFit1;
665 { efTRN, "-f", "eigenvec", ffREAD },
666 { efXVG, "-eig", "eigenval", ffOPTRD },
667 { efTPS, NULL, NULL, ffREAD },
668 { efNDX, NULL, NULL, ffOPTRD },
669 { efSTX, "-tar", "target", ffOPTRD},
670 { efSTX, "-ori", "origin", ffOPTRD},
671 { efEDI, "-o", "sam", ffWRITE }
673 #define NFILE asize(fnm)
674 edi_params.outfrq=100; edi_params.slope=0.0; edi_params.maxedsteps=0;
675 CopyRight(stderr,argv[0]);
676 parse_common_args(&argc,argv, 0 ,
677 NFILE,fnm,NPA,pa,asize(desc),desc,0,NULL,&oenv);
679 indexfile=ftp2fn_null(efNDX,NFILE,fnm);
680 EdiFile=ftp2fn(efEDI,NFILE,fnm);
681 TargetFile = opt2fn_null("-tar",NFILE,fnm);
682 OriginFile = opt2fn_null("-ori",NFILE,fnm);
685 for (ev_class=0; ev_class<evNr; ++ev_class) {
686 if (opt2parg_bSet(evOptions[ev_class],NPA,pa)) {
687 /*get list of eigenvectors*/
688 nvecs=sscan_list(&(listen[ev_class]),opt2parg_str(evOptions[ev_class],NPA,pa),evOptions[ev_class]);
689 if (ev_class<evStepNr-2) {
690 /*if apropriate get list of stepsizes for these eigenvectors*/
691 if (opt2parg_bSet(evStepOptions[ev_class],NPA,pa))
692 evStepList[ev_class]=
693 scan_vecparams(opt2parg_str(evStepOptions[ev_class],NPA,pa),evStepOptions[ev_class],nvecs);
694 else { /*if list is not given fill with zeros */
695 snew(evStepList[ev_class],nvecs);
696 for (i=0; i<nvecs; i++)
697 evStepList[ev_class][i]=0.0;
699 } else if (ev_class == evRADFIX && opt2parg_bSet(evStepOptions[ev_class],NPA,pa)) {
700 snew(evStepList[ev_class],nvecs);
701 for (i=0; i<nvecs; i++)
702 evStepList[ev_class][i]=radfix;
703 } else if (ev_class == evFLOOD) {
704 snew(evStepList[ev_class],nvecs);
706 /* Are we doing constant force flooding? In that case, we read in
707 * the fproj values from the command line */
708 if (opt2parg_bSet("-constF", NPA, pa))
710 evStepList[ev_class] = scan_vecparams(opt2parg_str("-constF",NPA,pa),"-constF",nvecs);
712 } else {}; /*to avoid ambiguity */
713 } else { /* if there are no eigenvectors for this option set list to zero */
714 listen[ev_class]=NULL;
715 snew(listen[ev_class],1);
716 listen[ev_class][0]=0;
720 /* print the interpreted list of eigenvectors - to give some feedback*/
721 for (ev_class=0; ev_class<evNr; ++ev_class) {
722 printf("Eigenvector list %7s consists of the indices: ",evOptions[ev_class]);
724 while(listen[ev_class][i])
725 printf("%d ",listen[ev_class][i++]);
730 EigvecFile=opt2fn("-f",NFILE,fnm);
732 /*read eigenvectors from eigvec.trr*/
733 read_eigenvectors(EigvecFile,&nav,&bFit1,
734 &xref1,&edi_params.fitmas,&xav1,&edi_params.pcamas,&nvec1,&eignr1,&eigvec1,&eigval1);
736 bTop=read_tps_conf(ftp2fn(efTPS,NFILE,fnm),
737 title,&top,&ePBC,&xtop,NULL,topbox,0);
741 printf("\nSelect an index group of %d elements that corresponds to the eigenvectors\n",nav);
742 get_index(atoms,indexfile,1,&i,&index,&grpname); /*if indexfile != NULL parameter 'atoms' is ignored */
744 gmx_fatal(FARGS,"you selected a group with %d elements instead of %d",
754 /* if g_covar used different coordinate groups to fit and to do the PCA */
755 printf("\nNote: the structure in %s should be the same\n"
756 " as the one used for the fit in g_covar\n",ftp2fn(efTPS,NFILE,fnm));
757 printf("\nSelect the index group that was used for the least squares fit in g_covar\n");
761 printf("\nNote: Apparently no fitting was done in g_covar.\n"
762 " However, you need to select a reference group for fitting in mdrun\n");
764 get_index(atoms,indexfile,1,&nfit,&ifit,&grpname);
767 copy_rvec(xtop[ifit[i]],xref1[i]);
775 if (opt2parg_bSet("-constF", NPA, pa))
777 /* Constant force flooding is special: Most of the normal flooding
778 * options are not needed. */
779 edi_params.flood.bConstForce = TRUE;
783 /* For normal flooding read eigenvalues and store them in evSteplist[evFLOOD] */
785 if ( listen[evFLOOD][0]!=0)
786 read_eigenvalues(listen[evFLOOD],opt2fn("-eig",NFILE,fnm),evStepList[evFLOOD],bHesse,kB*T);
788 edi_params.flood.tau=tau;
789 edi_params.flood.deltaF0=deltaF0;
790 edi_params.flood.deltaF=deltaF;
791 edi_params.presteps=eqSteps;
792 edi_params.flood.kT=kB*T;
793 edi_params.flood.bHarmonic=bHarmonic;
796 /* Trick: invert sign of Efl and alpha2 then this will give the same sign in the exponential and inverted sign outside */
797 edi_params.flood.constEfl=-constEfl;
798 edi_params.flood.alpha2=-sqr(alpha);
802 edi_params.flood.constEfl=constEfl;
803 edi_params.flood.alpha2=sqr(alpha);
809 /*number of system atoms */
810 edi_params.nini=atoms->nr;
813 /*store reference and average structure in edi_params*/
814 make_t_edx(&edi_params.sref,nfit,xref1,ifit );
815 make_t_edx(&edi_params.sav ,nav ,xav1 ,index);
818 /* Store target positions in edi_params */
819 if (opt2bSet("-tar",NFILE,fnm))
821 if (0 != listen[evFLOOD][0])
823 fprintf(stderr, "\nNote: Providing a TARGET structure has no effect when using flooding.\n"
824 " You may want to use -ori to define the flooding potential center.\n\n");
826 get_structure(atoms,indexfile,TargetFile,&edi_params.star,nfit,ifit,nav,index);
830 make_t_edx(&edi_params.star,0,NULL,index);
833 /* Store origin positions */
834 if (opt2bSet("-ori",NFILE,fnm))
836 get_structure(atoms,indexfile,OriginFile,&edi_params.sori,nfit,ifit,nav,index);
840 make_t_edx(&edi_params.sori,0,NULL,index);
844 write_the_whole_thing(ffopen(EdiFile,"w"), &edi_params, eigvec1,nvec1, listen, evStepList);