NAMD
SimParameters.h
Go to the documentation of this file.
1 
7 /*****************************************************************************
8  * $Source: /home/cvs/namd/cvsroot/namd2/src/SimParameters.h,v $
9  * $Author: jim $
10  * $Date: 2017/03/30 20:06:17 $
11  * $Revision: 1.1248 $
12  *****************************************************************************/
13 
14 #ifndef SIMPARAMETERS_H
15 #define SIMPARAMETERS_H
16 
17 #include "common.h"
18 #include "Vector.h"
19 #include "Lattice.h"
20 #include "imd.h"
21 
22 #include "MGridforceParams.h"
23 #include "GroupRestraintsParam.h"
24 
25 class ParseOptions;
26 class Communicate;
27 class ConfigList;
28 class MIStream;
29 class MOStream;
30 
31 // The class SimParameters is really just a glorified structure used to
32 // maintain the global simulation parameters. The only functions
33 // associated with the class are used to get the parameters from the
34 // ConfigList object, to send that Parameters from the master node
35 // to the other nodes, and to receive the Parameters on the other nodes.
36 
37 
38 // The following definitions are used to distinguish between possible
39 // bonded exclusion settings
40 typedef int ExclusionSettings;
41 
42 #define NONE 0
43 #define ONETWO 1
44 #define ONETHREE 2
45 #define ONEFOUR 3
46 #define SCALED14 4
47 
48 // The following definitions are used to distinguish between multiple
49 // timestep integration schemes
50 typedef int MTSChoices;
51 
52 #define NAIVE 0
53 #define VERLETI 1
54 
55 // The following definitions are used to distinuish between multiple
56 // long-short range force splittings
57 #define SHARP 0
58 #define XPLOR 1
59 #define C1 2
60 #define C2 3
61 
62 // The following definitions are used to distinguish among load
63 // balancers and their strategies
64 #define LDBAL_NONE 0
65 #define LDBAL_CENTRALIZED 1 // default
66 #define LDBAL_HYBRID 2
67 
68 #define LDBSTRAT_DEFAULT 10 // default
69 #define LDBSTRAT_COMPREHENSIVE 11
70 #define LDBSTRAT_REFINEONLY 12
71 #define LDBSTRAT_OLD 13
72 
73 // The following definitions are used to distinguish between patch-splitting
74 // strategies
75 #define SPLIT_PATCH_POSITION 0 // atom position determines patch
76 #define SPLIT_PATCH_HYDROGEN 1 // hydrogen groups are not broken up
77 
78 // The following definitions are used to distinguish the range of rigid
79 // bond calculations: none, all bonds to hydrogen, or only water
80 #define RIGID_NONE 0
81 #define RIGID_ALL 1
82 #define RIGID_WATER 2
83 
84 // Added by JLai -- The following definitions are used to distinguish
85 // the different GoMethodologies available to the Go program
86 // -- 6.3.11
87 typedef int GoChoices;
88 #define GO_MATRIX 1
89 #define GO_SPARSE 2
90 #define GO_LOWMEM 3
91 
92 // Used for controlling PME parallelization with ckloop
93 // The higher level will include all parallelization for lower ones
94 // E.g. If setting useCkLoop to 3, then xpencil's kspace, all
95 // backward ffts and send_untrans/ungrid routines will be parallelized
96 #define CKLOOP_CTRL_PME_UNGRIDCALC 6
97 #define CKLOOP_CTRL_PME_FORWARDFFT 5
98 #define CKLOOP_CTRL_PME_SENDTRANS 4
99 #define CKLOOP_CTRL_PME_KSPACE 3
100 #define CKLOOP_CTRL_PME_BACKWARDFFT 2
101 #define CKLOOP_CTRL_PME_SENDUNTRANS 1
102 
103 // These macros are used for routines of writing atom positions and velocities
104 // when NAMD crashes. To determine if we need to write the atoms to files, we
105 // can use code like the following example:
106 // if (simParams->crashOutputFlag & NAMD_CRASH_ATOM_TOO_FAST) {
107 // write_atoms_to_file();
108 // }
109 // Currently only the routine for errors like "atoms are moving too fast" is
110 // implemented, so there is only one macro. We may add more macros in the future
111 // to handle different type of errors.
112 #define NAMD_CRASH_ATOM_TOO_FAST 0x1
113 // In the future if more macros like NAMD_CRASH_XXX are added, then
114 // the follwing line should be changed to
115 // #define NAMD_CRASH_ALL (NAMD_CRASH_ATOM_TOO_FAST + NAMD_CRASH_XXX)
116 #define NAMD_CRASH_ALL NAMD_CRASH_ATOM_TOO_FAST
117 
118 // Bitmask for calculating bonded interactions on GPU.
119 #define NAMD_BONDEDGPU_BONDS (1 << 0)
120 #define NAMD_BONDEDGPU_ANGLES (1 << 1)
121 #define NAMD_BONDEDGPU_DIHEDRALS (1 << 2)
122 #define NAMD_BONDEDGPU_IMPROPERS (1 << 3)
123 #define NAMD_BONDEDGPU_EXCLS (1 << 4)
124 #define NAMD_BONDEDGPU_CROSSTERMS (1 << 5)
125 #define NAMD_BONDEDGPU_THOLES (1 << 6)
126 #define NAMD_BONDEDGPU_ANISOS (1 << 7)
127 #define NAMD_BONDEDGPU_ONEFOURENBTHOLES (1 << 8)
128 #define NAMD_BONDEDGPU_ALL ( \
129  NAMD_BONDEDGPU_BONDS +\
130  NAMD_BONDEDGPU_ANGLES +\
131  NAMD_BONDEDGPU_DIHEDRALS +\
132  NAMD_BONDEDGPU_IMPROPERS +\
133  NAMD_BONDEDGPU_EXCLS +\
134  NAMD_BONDEDGPU_CROSSTERMS +\
135  NAMD_BONDEDGPU_THOLES +\
136  NAMD_BONDEDGPU_ANISOS +\
137  NAMD_BONDEDGPU_ONEFOURENBTHOLES)
138 
140 {
141 private:
142 public:
143 
144 // MAKE SURE THAT THIS CLASS CAN BE BIT COPIED OR YOU WILL HAVE TO
145 // ADD SPECIAL CODE TO send_SimParameters() and receive_SimParameters()
146 
147  int mshakeOn;
148  int lincsOn;
149 
150 #if defined(NAMD_NVTX_ENABLED) || defined(NAMD_CMK_TRACE_ENABLED) || defined(NAMD_ROCTX_ENABLED)
151  int beginEventPatchID;
152  int endEventPatchID;
153  int beginEventStep;
154  int endEventStep;
155 #endif
156 
157 #ifdef TIMER_COLLECTION
158  double timerBinWidth; // default 1
159 #endif
160 
166  Bool SOAintegrateOn; // use SOA integration routine for higher performance
167 
179  // this is set true to indicate dynamics is being performed
180 
182  // window size for moving averages of reductions for GPU-resident mode
183 
184  Bool nsPerDayOn; // prints ns/day instead of days/ns
185 
186 
187  Bool lonepairs; // enable lone pairs
188  WaterModel watmodel; // integer code for the water model in use
189  // choices are defined in common.h
190  Bool LJcorrection; // flag for whether water tail corrections should be used
191  Bool LJcorrectionAlt; // flag for whether alternative tail corrections should be used
192  BigReal dt; // Timestep size
193  int N; // Number of steps to be performed
194  int stepsPerCycle; // Number of timesteps per cycle
195 
196  zVector cellBasisVector1; // Basis vector for periodic cell
197  zVector cellBasisVector2; // Basis vector for periodic cell
198  zVector cellBasisVector3; // Basis vector for periodic cell
199  zVector cellOrigin; // Fixed center of periodic cell
200  Lattice lattice; // All data for periodic cell
201 
202  int nonbondedFrequency; // Number of timesteps between
203  // nonbonded evaluation
204  int fullElectFrequency; // Number of timesteps between
205  // full electrostatic evaluation
206  int fullDispersionFrequency; // Number of timesteps between
207  // full LJ dispersion evaluation
208  BigReal fmaTheta; // DPMTA theta value
209  int ldBalancer; // None, Centralized or Hybrid
210  int ldbStrategy; // What load balancing strategy to use
211  int ldbPeriod; // How often to do load balancing
212  int firstLdbStep; // What step to do the first
213  // load-balance on.
214  int lastLdbStep; // What step to do the last
215  // load-balance on.
216  int hybridGroupSize; // hybrid group size
217  BigReal ldbBackgroundScaling; // scaling factor for background load
218  BigReal ldbPMEBackgroundScaling;// scaling factor for PME background
219  BigReal ldbHomeBackgroundScaling;// scaling factor for home background
220  BigReal ldbRelativeGrainsize; // fraction of average load per compute
221 
222  int traceStartStep; //the timestep when trace is turned on, default to 3*firstLdbStep;
223  int numTraceSteps; //the number of timesteps that are traced, default to 2*ldbPeriod;
224 
225 #ifdef MEASURE_NAMD_WITH_PAPI
226  Bool papiMeasure; //default to false
227  int papiMeasureStartStep; //the timestep when to measure using PAPI, default to 3*firstLdbStep;
228  int numPapiMeasureSteps; //the number of timesteps when performance are measured with PAPI, default to 40;
229 #endif
230 
231  Bool outputMaps; //control whether to dump compute/patch map before load balancing
232  Bool simulateInitialMapping; //if true, the initial mapping during startup is dumped and exit
235  Bool disableTopology; // ignore torus information during patch placement
236  Bool verboseTopology; // print torus information during patch placement
237 
238  Bool benchTimestep; //only cares about benchmarking the timestep, so no file output to save SUs for large-scale benchmarking
239 
240  //whether to use CkLoop library to parallelize a loop in a function like OpenMP.
241  //It has multiple control levels. The higher the value is (must be positive), the more parallelization will be performed
242  //Currently, it is mainly used for PME computation. The default value is 0, meaning it is disabled
243  //Refer to macros CKLOOP_CTRL_* in this file for the ordering of different levels
244  int useCkLoop;
245 
246  int twoAwayX; // half-size patches in X dimension
247  int twoAwayY; // half-size patches in Y dimension
248  int twoAwayZ; // half-size patches in Z dimension
249  int maxPatches; // maximum patch count
250  Bool ldbUnloadPME; // unload processors doing PME
251  Bool ldbUnloadZero; // unload processor 0
252  Bool ldbUnloadOne; // unload processor 1
253  Bool ldbUnloadOutputPEs; // unload output processors
254  Bool noPatchesOnZero; // no patches on processor 0
255  Bool noPatchesOnOutputPEs; // no patches on output PEs
256  Bool noPatchesOnOne; // no patches on processor 1
257 
258  BigReal initialTemp; // Initial temperature for the
259  // simulation
260  Bool comMove; // Should the center of mass be
261  // able to move
262  Bool zeroMomentum; // remove momentum drift from PME
263  Bool zeroMomentumAlt; // alternate method for testing
264  Bool wrapWater; // Wrap water around on output
265  Bool wrapAll; // Wrap clusters around on output
266  Bool wrapNearest; // Wrap to closest image to origin
267  BigReal dielectric; // Dielectric constant
268  ExclusionSettings exclude; // What electrostatic exclusions should
269  // be made
270 
271  // scale14alt will override if scale14 is not set
272  BigReal scale14; // Scaling factor "1-4scaling"
273  // for 1-4 electrostatics
274  BigReal scale14alt; // Alternatively named sim parameter
275  // "oneFourScaling" so that this can
276  // be set using scripting language
277 
278  BigReal nonbondedScaling; // Scaling factor for nonbonded forces
279  int dcdFrequency; // How often (in timesteps) should
280  // a DCD trajectory file be updated
281  int dcdUnitCell; // Whether to write unit cell information in the DCD
282  int velDcdFrequency; // How often (in timesteps) should
283  // a velocity DCD file be updated
284  int forceDcdFrequency; // How often (in timesteps) should
285  // a force DCD file be updated
286  int xstFrequency; // How often (in timesteps) should
287  // a XST trajectory file be updated
288  int dcdSelectionOn; // enable dcd selection options
289  char auxFilename[NAMD_FILENAME_BUFFER_SIZE]; // auxilary output filename
290  char dcdFilename[NAMD_FILENAME_BUFFER_SIZE]; // DCD filename
291  char velDcdFilename[NAMD_FILENAME_BUFFER_SIZE]; // Velocity DCD filename
292  char forceDcdFilename[NAMD_FILENAME_BUFFER_SIZE]; // Force DCD filename
293  char xstFilename[NAMD_FILENAME_BUFFER_SIZE]; // Extended system trajectory filename
294  char outputFilename[NAMD_FILENAME_BUFFER_SIZE]; // Output file name. This name will
295  // have .coor appended to it
296  // for the coordinates and
297  // .vel appended to
298  // it for the velocities
299  char restartFilename[NAMD_FILENAME_BUFFER_SIZE]; // Base name of the restart file
300  char crashFilename[NAMD_FILENAME_BUFFER_SIZE]; // Base name of the crash CSV file containing atom positions and velocities
301  int crashOutputFlag; // Bitfield to control the crash output
302  // std::vector <DCDParams> dcdUserDefined; //parameters for user defined DCD lists
303  int restartFrequency; // How often (in timesteps) shoud the
304  // restart files be updated
305  Bool restartSave; // unique filenames for restart files
306  Bool restartSaveDcd; // unique filenames for DCD files
307  Bool binaryRestart; // should restart files be
308  // binary format rather than PDB
309  Bool binaryOutput; // should output files be
310  // binary format rather than PDB
311  BigReal cutoff; // Cutoff distance
312  BigReal margin; // Fudge factor on patch size
313  BigReal patchDimension; // Dimension of each side of a patch
314  // This is either cutoff+margin or
315  // pairlistDist+margin depending on
316  // whether or not switching is on
317  // or not
318  BigReal limitDist; // Distance below which nonbonded
319  // forces between atoms are limited
320  Bool switchingActive; // Flag TRUE->using switching function
321  // for electrostatics and vdw
322  Bool vdwForceSwitching; // Flag TRUE->using force switching
323  // function for vdw
324  BigReal switchingDist; // Distance at which switching
325  // becomes active
326  Bool martiniSwitching; // Flag TRUE->use Martini residue-based
327  // coarse-grain switching function
328  Bool martiniDielAllow; // Allow non-standard dielectric constant
329  // for use with Martini when dielectric != 15.0
330  BigReal pairlistDist; // Distance within which atom pairs
331  // should be added to pairlist
332  int pairlistMinProcs; // Minimum number of processors
333  // to enable pairlists
334  int usePairlists; // Derived from pairlistMinProcs
335 
336  int pairlistsPerCycle; // regenerate x times per cycle
337  BigReal pairlistShrink; // tol *= (1 - x) on regeneration
338  BigReal pairlistGrow; // tol *= (1 + x) on trigger
339  BigReal pairlistTrigger; // trigger is atom > (1 - x) * tol
340  int outputPairlists; // print pairlist warnings this often
341 
342  Bool constraintsOn; // Flag TRUE-> harmonic constraints
343  // active
344  int constraintExp; // Exponent for harmonic constraints
345 
346  /* BEGIN gf */
347  Bool gridforceOn; // Flag TRUE -> gridforce active
348  Bool gridforceVolts; // Flag TRUE -> gridforce using volts as units
349  zVector gridforceScale; // Gridforce scale factor
350  Bool gridforceContA1; // Flag TRUE -> grid continuous in A1 direction
351  Bool gridforceContA2; // Flag TRUE -> grid continuous in A2 direction
352  Bool gridforceContA3; // Flag TRUE -> grid continuous in A3 direction
353  zVector gridforceVOffset; // Gridforce potential offsets
354  Bool gridforceLite; // Flag TRUE -> use lightweight, fast, feature-poor gridforce
355  Bool gridforcechecksize; //Flag TRUE -> check if grid is larger than PBC cell dimensions
356  /* END gf */
359 
360  /* Begin Group restraints */
361  Bool groupRestraintsOn; // Turn on the group restraints
362  GroupRestraintList groupRestraints; // To store restraints parameters for each group
363  int groupRestraintsCount; // count of how many groups
364  /* End Group restraints */
365 
366  //****** BEGIN selective restraints (X,Y,Z) changes
367  Bool selectConstraintsOn; // Flag TRUE-> selective restraints
368  // active
370  constrZOn; // Flag TRUE-> select which Cartesian
371  // component to restrain
372  //****** END selective restraints (X,Y,Z) changes
373 
374  // spherical constraints
377 
378  BigReal constraintScaling; // Scaling factor for constraint forces
379 
380  //****** BEGIN CHARMM/XPLOR type changes
381  Bool paraTypeXplorOn; // FLAG TRUE-> parametrs are XPLOR format (default)
382  Bool paraTypeCharmmOn; // FLAG TRUE-> parametrs are CHARMM format
383  //****** END CHARMM/XPLOR type changes
384 
385  // Ported by JLai -- JE - Go
386  Bool goGroPair; // FLAG FALSE->Explicit Gromacs pairs will be calculated
387  Bool goForcesOn; // FLAG TRUE-> Go forces will be calculated
388  char goParameters[NAMD_FILENAME_BUFFER_SIZE]; // File for Go parameters
389  char goCoordinates[NAMD_FILENAME_BUFFER_SIZE]; // File for Go structure and atom chain types
390  //JLai 6.3.11
391  GoChoices goMethod; // Integer for Go method -- 1) Matrix-Go, 3) Low-mem-Go
392  // End of port -- JL
393 
394  //****** BEGIN moving constraints changes
395  Bool movingConstraintsOn; // Flag TRUE-> moving constraints
396  // active
397  zVector movingConsVel; // Velocity of the movement, A/timestep
398  //****** END moving constraints changes
399  //****** BEGIN rotating constraints changes
400  Bool rotConstraintsOn; // Flag TRUE-> rotating constraints
401  // active
402  zVector rotConsAxis; // Axis of rotation
403  zVector rotConsPivot; // Pivot point of rotation
404  BigReal rotConsVel; // Velocity of rotation, Deg/timestep
405  //****** END rotating constraints changes
406 
407  //****** BEGIN moving drag changes
408  Bool movDragOn; // Flag TRUE-> moving drag active
409  char movDragFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file defining dragged atoms
410  // by non-zero value in the column
411  BigReal movDragGlobVel; // global drag velocity (A/step)
412  char movDragVelFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ scale moving drag
413  // velocity for each atom
414  //****** END moving drag changes
415  //****** BEGIN rotating drag changes
416  Bool rotDragOn; // Flag TRUE-> rotating drag active
417  char rotDragFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file defining dragged atoms
418  // by non-zero value in the column
419  char rotDragAxisFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ define axes for atoms;
420  char rotDragPivotFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ define pivots for atoms
421  BigReal rotDragGlobVel; // global drag velocity (deg/step)
422  char rotDragVelFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; B or O scales angular
423  // velocity for each atom
424  //****** END rotating drag changes
425  //****** BEGIN "constant" torque changes
426  Bool consTorqueOn; // Flag TRUE-> "constant" torque active
427  char consTorqueFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file defining torqued atoms
428  // by non-zero value in the column
429  char consTorqueAxisFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ define axes for atoms;
430  char consTorquePivotFile[NAMD_FILENAME_BUFFER_SIZE];// PDB file; XYZ define pivots for atoms
431  BigReal consTorqueGlobVal; // global "torque" (Kcal/(mol*A^2))
432  char consTorqueValFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; B or O scales "torque"
433  // for each atom
434  //****** END "constant" torque changes
435 
436  //****** BEGIN SMD constraints changes
437  Bool SMDOn; // Flag TRUE-> SMD constraints active
438  BigReal SMDVel; // Velocity of the movement, A/timestep
439  zVector SMDDir; // Direction of the movement
440  BigReal SMDk; // Elastic constant for SMD
441  BigReal SMDk2; // Transverse elastic constant for SMD
442  char SMDFile[NAMD_FILENAME_BUFFER_SIZE]; // File for SMD information
443  int SMDOutputFreq; // Output frequency for SMD constr.
444  //****** END SMD constraints changes
445 
446  //****** BEGIN tabulated energy section
450  char tableInterpType[128];
453  //****** END tabulated energy section
454 
455  // TMD
462 
463  //Symmetry restraints
470 
471 
472 //fepb
473  Bool alchOnAtStartup; // Ensure that alchemy is set up properly
476  Bool alchOn; // Doing alchemical simulation?
477  Bool alchFepOn; // Doing alchemical simulation?
478  Bool singleTopology; // Using single topology setup?
479  Bool sdScaling; // Scaling S-D bond terms in single topology?
480  Bool alchThermIntOn; // Doing thermodynamic integration?
481  Bool alchWCAOn; // Using WCA decomposition for vdWs?
482  int alchMethod; // Which alchemical method to use? fep or ti
483  BigReal alchLambda; // lambda for dynamics
484  BigReal alchLambda2; // lambda for comparison
485  BigReal alchLambdaIDWS; // alternate lambda for interleaved double-wide sampling
486  int alchIDWSFreq; // freq with which lambda2 changes to lambdaIDWS
487  int alchLambdaFreq; // freq. (in steps) with which lambda changes
488  // from alchLambda to alchLambda2
489  BigReal getCurrentLambda(const int) const; // getter for changing lambda
490  BigReal getCurrentLambda2(const int) const; // getter for alternating lambda2 in IDWS
491  int setupIDWS(); // activates IDWS and sets alchIDWSFreq
492  BigReal getLambdaDelta(void) const; // getter for lambda increment
493  BigReal alchTemp; // temperature for alchemical calculation
494  int alchOutFreq; // freq. of alchemical output
495  Bool alchEnsembleAvg; //if do ensemble average for the net free energy difference
496  char alchOutFile[NAMD_FILENAME_BUFFER_SIZE]; // alchemical output filename
497  int alchEquilSteps; // # of equil. steps in the window
498  BigReal alchVdwShiftCoeff; // r2 shift coeff used for generating
499  // the alchemical altered vdW interactions
500  BigReal alchElecLambdaStart; // lambda value for starting point of
501  // electrostatic interactions of
502  // exnihilated particles. For annihilated
503  // particles the starting point is
504  // (1-alchElecLambdaStart)
505  BigReal getElecLambda(const BigReal) const; // return min[0,x/(1-elecStart)]
506  BigReal alchVdwLambdaEnd; // lambda value for endpoint of vdW
507  // interactions of exnihilated particles.
508  // For annihilated particles the endpoint is
509  // (1-alchVdwLambdaEnd)
510  BigReal getVdwLambda(const BigReal) const; // return max[1,x/vdwEnd]
511  BigReal alchRepLambdaEnd; // lambda value for endpoint of repulsive vdW
512  // interactions of exnihilated particles.
513  // For annihilated particles the endpoint is
514  // (1-alchRepLambdaEnd). This also implies the
515  // START for attractive vdW interactions.
516  BigReal getRepLambda(const BigReal) const; // return max[1,x/repEnd]
517  BigReal alchBondLambdaEnd; // lambda value for endpoint of bonded
518  // interactions involving exnihilated particles.
519  // For annihilated particles the endpoint is
520  // (1-alchBondLambdaEnd)
521  BigReal getBondLambda(const BigReal) const; // return max[1,x/bondEnd]
522  Bool alchDecouple; // alchemical decoupling rather than annihilation
523  Bool alchBondDecouple; // decouple purely alchemical bonds
524  size_t alchGetNumOfPMEGrids() const; // get the number of PME grids required by alchemy
525 //fepe
526 
527 
528  Bool lesOn; // Locally enhanced sampling?
529  int lesFactor; // local enhancement factor
530  Bool lesReduceTemp; // Reduce enhanced atom temperature?
531  Bool lesReduceMass; // Reduce enhanced atom mass?
532 
533  // REST2
547  Bool extForcesOn; // Are ext command forces present?
551 
552 
553  // Defines variables for QM/MM calculations
554  Bool qmForcesOn; // Are QM/MM command forces present?
559  char qmSoftware[128];
560  char qmChrgModeS[16];
562  char qmColumn[16];
568  int qmFormat ;
571 
573  char qmBondColumn[16];
579  char qmBondSchemeS[16] ;
582  char qmPCSwitchTypeS[16];
584  char qmPCSchemeS[16];
587 
593 
595  int qmLSSFreq ;
596  char qmLSSResname[5] ;
597  char qmLSSModeS[16];
599 
602 
604  int qmOutFreq ;
606 
607  Bool printBadContacts; //print indices of bad contacts being moved downhill
608 
609  //gbis implicit solvent parameters
610  Bool GBISOn; //do generalized born implicit solvent
612  Bool GBISserOn; //do generalized born implicit solvent serial
615  BigReal kappa; //debye screening length; k = sqrt(ion concentration mol/L ) / 0.304
617  BigReal gbis_delta; //three parameters for born radius calc
620  BigReal alpha_cutoff; //pairwise cutoff for integrating born radius
621  BigReal alpha_max; //maximum allowable born radius
622  Bool LCPOOn; //do LCPO SASA for GBSA
623  BigReal surface_tension; //surface tension (kcal/mol/Ang^2) for LCPO
624 
625  Bool drudeOn; // Perform integration of Drude oscillators?
626  Bool drudeHardWallOn; // Apply maximum Drude bond length restriction?
627  BigReal drudeTemp; // (low) temperature for freezing Drude oscillators
628  BigReal drudeDamping; // Langevin damping coefficient (1/ps)
629  // defaults to langevinDamping
630  BigReal drudeBondLen; // Length beyond which to apply quartic
631  // restraining potential to Drude bond
632  BigReal drudeBondConst; // Force constant for restraining potential
633  BigReal drudeNbtholeCut; // Radius of thole pair interaction
634 
635  Bool pairInteractionOn; // Calculate pair interactions?
636  int pairInteractionGroup1; // Interaction group 1.
637  int pairInteractionGroup2; // Interaction group 2.
638  Bool pairInteractionSelf; // Compute just within group.
639 
640  Bool cosAngles; // Can some angles be cos-based
641  Bool globalForcesOn; // Are global forces present?
642  Bool tclForcesOn; // Are Tcl forces present?
650 #ifdef NAMD_TCL
651  Bool tclIsThreaded; // Is Tcl library thread-safe?
652 #endif
653  Bool tclBCOn; // Are Tcl boundary forces present
654  char *tclBCScript; // Script defining tclBC calcforces
655  char tclBCArgs[NAMD_FILENAME_BUFFER_SIZE]; // Extra args for calcforces command
656  Bool freeEnergyOn; // Doing free energy perturbation?
657  Bool miscForcesOn; // Using misc forces?
658  Bool colvarsOn; // Using the colvars module?
659 
660  Bool fixedAtomsOn; // Are there fixed atoms?
661  Bool fixedAtomsForces; // Calculate forces anyway?
662  Bool fixedAtomsForceOutput; // Output fixed forces?
663 
664  Bool langevinOnAtStartup; // Ensure that langevin is set up properly
665  Bool langevinOn; // Flag TRUE-> langevin dynamics active
666  BigReal langevinTemp; // Temperature for Langevin dynamics
667  BigReal langevinDamping; // Damping coefficient (1/ps)
668  Bool langevinHydrogen; // Flag TRUE-> apply to hydrogens
673  Bool langevin_useBAOAB; // Flag TRUE-> use the experimental BAOAB integrator for NVT instead of the BBK one
674  // See Leimkuhler and Matthews (AMRX 2012); implemented in NAMD by CM June2012
675 
676  // BEGIN LA
677  Bool loweAndersenOn; // Flag TRUE-> Lowe-Andersen dynamics active
678  BigReal loweAndersenTemp; // Temperature for Lowe-Andersen dynamics
679  BigReal loweAndersenRate; // Collision frequency for Lowe-Andersen dynamics (1/ps)
680  BigReal loweAndersenCutoff; // Cutoff radius for Lowe-Andersen dynamics
681  // END LA
682 
683  Bool globalOn; // Flag TRUE-> use global integrator
684  Bool dihedralOn; // Flag TRUE-> dihedral dynamics active
685  Bool COLDOn; // Flag TRUE-> constrained overdamped
686  // langevin dynamics active
687  BigReal COLDRate; // Damping coefficient for COLD.
688  BigReal COLDTemp; // Temperature for COLD.
689 
690  Bool tCoupleOn; // Flag TRUE-> Temperature coupling
691  // active
692  BigReal tCoupleTemp; // Temperature for temp coupling
693 
715  int rescaleFreq; // Velocity rescale frequency
716  BigReal rescaleTemp; // Temperature to rescale to
717 
718  Bool accelMDOn; // Perform accelerated MD
719  Bool accelMDdihe; // Apply boost to the dihedral potential
720  Bool accelMDdual; // dual boost mode
721  Bool accelMDDebugOn; // Debugging accelerated MD
722  BigReal accelMDFirstStep; // First aMD step
723  BigReal accelMDLastStep; // Last aMD step
724  int accelMDOutFreq; // aMD output frequency
725  BigReal accelMDE; // aMD E
726  BigReal accelMDalpha; // aMD alpha
727  BigReal accelMDTE; // E for total potential in the dual boost mode
728  BigReal accelMDTalpha; // alpha for total potential in the dual boost mode
729 
730  Bool accelMDG; // Perform Gaussian accelMD calculation
731  int accelMDGiE; // Flag to set the mode iE in Gaussian accelMD
732  int accelMDGcMDSteps; // Number of cMD steps
733  int accelMDGEquiSteps; // Number of quilibration steps after adding boost potential
734  int accelMDGcMDPrepSteps; // Number of preparation cMD steps
735  int accelMDGEquiPrepSteps; // Number of preparation equilibration steps
736  int accelMDGStatWindow; // Number of steps to calc avg and std
737  BigReal accelMDGSigma0P; // upper limit of std of total potential
738  BigReal accelMDGSigma0D; // upper limit of std of dihedral potential
739  Bool accelMDGRestart; // Flag to set use restart file in Gaussian accelMD
740  char accelMDGRestartFile[NAMD_FILENAME_BUFFER_SIZE]; // restart file name
741  Bool accelMDGresetVaftercmd; // Flag to reset potential after first accelMDGcMDSteps steps
742 
743  /* Begin Adaptive Temperature Sampling */
744  Bool adaptTempOn; // is adaptTempOn
745  Bool adaptTempDebug; // Debuggin adaptive temperature sampling
746  int adaptTempFirstStep; // First adaptTemp step
747  int adaptTempLastStep; // Last adaptTemp step
748  int adaptTempOutFreq; // adaptTemp output frequency
749  int adaptTempFreq; // Steps between adaptTemp updates
750  BigReal adaptTempTmin; // Lower temperature bound
751  BigReal adaptTempTmax; // Upper temperature bound
752  BigReal adaptTempAutoDt; // Auto jump size. Value determines upper bound, adaotTempDt determines lower bound
753  int adaptTempBins; // Number of bins to store average energy values
754  BigReal adaptTempDt; // timestep for adaptTemp updates - only affects Temperature random walk
755  BigReal adaptTempCgamma; // Cgamma variable for adaptive bin averaging Cgamma = 0 is normal Averaging. 1 > Cgamma >= 0
756  Bool adaptTempLangevin; // Couple to Langevin Thermostat
757  Bool adaptTempRescale; // Couple to Vel. Rescaling
758  char adaptTempInFile[NAMD_FILENAME_BUFFER_SIZE]; // Restart information for adaptTemp to read
759  char adaptTempRestartFile[NAMD_FILENAME_BUFFER_SIZE]; // File to write restart information
760  int adaptTempRestartFreq; // Frequency of writing restart output
761  Bool adaptTempRandom; // Do we assign random temperatures when we step out of [Tmin,Tmax]?
762  /* End Adaptive Temperature Sampling */
763 
764  int reassignFreq; // Velocity reassignment frequency
765  BigReal reassignTemp; // Temperature to reassign to
766  BigReal reassignIncr; // Added to reassignTemp each time
767  BigReal reassignHold; // Hold reassignTemp at this value
768 
769  Bool useGroupPressure; // Use group rather than atomic
770  // quantities for pressure calc
771 
772  Bool excludeFromPressure; // Flag TRUE-> some atoms not rescaled
773 
774  Bool useFlexibleCell; // Use anisotropic cell fluctuations
775  Bool useConstantArea; // x,y dimensions fixed.
776  Bool useConstantRatio; // x,y ratio fixed.
777 
778  Bool fixCellDims; // fix the cell dimensions
782 
783  Bool berendsenPressureOn; // Berendsen pressure bath
788 
789  Bool langevinPistonOn; // Langevin piston pressure control
790  Bool langevinPistonBarrier; // Turn off to extrapolate cell
795 
796  Bool monteCarloPressureOnAtStartup; // Ensure that MC pressure is set up properly
797  Bool monteCarloPressureOn; // MonteCarlo pressure control
798  BigReal monteCarloPressureTarget; // target pressure
799  BigReal monteCarloTemp; // Temperature used in Monte Carlo acceptance critiria
800  BigReal monteCarloAcceptanceRate; // target acceptance rate
801  zVector monteCarloMaxVolume; // initial maximum volume change for each dimension
802  int monteCarloAdjustmentFreq; // frequency of adjusting maximum scaling factor
803  int monteCarloPressureFreq; // frequency of applying MC pressure
804 
805  Bool multigratorOn; // Multigrator temperature and/or pressure control
813 
815 
816  Bool pressureProfileOn; // Compute lateral pressure profile?
817  int pressureProfileSlabs; // Number of slabs
818  int pressureProfileFreq; // How often to store profile data
820  Bool pressureProfileEwaldOn; // Compute Ewald contribution?
824 
826  zVector strainRate2; // off diagonal elements (xy, xz, yz)
827 
828  unsigned int randomSeed; // Seed for random number generator
829 
830  Bool FMAOn; // Flag TRUE-> FMA active
831  int FMALevels; // Number of Levels for FMA
832  int FMAMp; // Number of multipole terms for FMA
833  Bool FMAFFTOn; // FFT on/off flag for FMA
834  int FMAFFTBlock; // FFT blocking factor for FMA
835 
836  Bool fullDirectOn; // Should direct calculations of
837  // full electrostatics be performed?
838 
839  Bool MSMOn; // enable MSM (multilevel summation method)
840  // for long-range electrostatics
841 
842  int MSMQuality; // choose MSM quality 0 (low) - 3 (high), using
843  // optimal combination of approximation and splitting
844  // defaults to "low" for fastest performance
845 
846  int MSMApprox; // choose MSM approximation
847  // defaults to "cubic" (low) for fastest performance
848 
849  int MSMSplit; // choose MSM splitting function
850  // defaults to "Taylor2" (low) for fastest performance
851 
852  int MSMLevels; // select number of MSM levels
853  // default (0) adapts number of levels to the
854  // system for fastest performance
855 
856  int MSMBlockSizeX; // controls size of parallel work decomposition
857  int MSMBlockSizeY; // controls size of parallel work decomposition
858  int MSMBlockSizeZ; // controls size of parallel work decomposition
859 
860  BigReal MSMGridSpacing; // defaults to 2.5 A, best for atomic systems
861 
862  BigReal MSMPadding; // pad grid along non-periodic boundaries
863  // defaults to 2.5 A
864  // increase if atoms are drifting beyond
865  // edge of grid, which will terminate
866  // simulation prematurely
867 
868  BigReal MSMxmin; // define extent of non-periodic boundaries
874 
875  Bool MsmSerialOn; // use serial MSM solver for testing
876 
880 
881  // LJ-PME parameters
882  Bool LJPMEOn; // Flag TRUE -> LJ-PME active
883  BigReal LJPMETolerance; // LJ Direct space tolerance
884  BigReal LJPMEEwaldCoefficient; // From tolerance and cutoff
885  int LJPMEInterpOrder; // LJ-PME Order of interpolation
886  int LJPMEGridSizeX; // No. of grid points in x dim
887  int LJPMEGridSizeY; // No. of grid points in y dim
888  int LJPMEGridSizeZ; // No. of grid points in z dim
889  BigReal LJPMEGridSpacing; // Maximum spacing between points
890  //Bool LJFFTWUseWisdom; // Read/save wisdom file for FFTW
891  //char LJFFTWWisdomFile[NAMD_FILENAME_BUFFER_SIZE]; // Wisdom file name for FFTW
892  //char *LJFFTWWisdomString;
893  Bool LJPMESerialOn; // Use serial implementation only
894  Bool LJPMESerialRealSpaceOn; // also use (slow) real space part
895 
896  // PME Parameters
897  Bool PMEOn; // Flag TRUE -> PME active
898  BigReal PMETolerance; // Direct space tolerance
899  BigReal PMEEwaldCoefficient; // From tolerance and cutoff
900  int PMEInterpOrder; // Order of interpolation
901  int PMEGridSizeX; // No. of grid points in x dim
902  int PMEGridSizeY; // No. of grid points in y dim
903  int PMEGridSizeZ; // No. of grid points in z dim
904  BigReal PMEGridSpacing; // Maximum spacing between points
905  int PMEProcessors; // No. of processors to use
906  int PMEMinSlices; // Min slices per PME slab
907  int PMEMinPoints; // Min points per PME pencil
908  Bool PMEBarrier; // Use barrier before sendTrans
909  int PMEPencils; // Size of pencil grid in each dim
910  int PMEPencilsX; // Size of pencil grid in X dim
911  int PMEPencilsY; // Size of pencil grid in Y dim
912  int PMEPencilsZ; // Size of pencil grid in Z dim
913  int PMEPencilsYLayout; // Y pencil layout strategy
914  int PMEPencilsXLayout; // X pencil layout strategy
915  int PMESendOrder; // Message ordering strategy
916  Bool PMEOffload; // Offload reciprocal sum to accelerator
917 
918  Bool useDPME; // Flag TRUE -> old DPME code
919 
948 
949  #ifdef OPENATOM_VERSION
950  Bool openatom; // Flag TRUE -> OpenAtom QM/MM active
951  #endif // OPENATOM_VERSION
952 
953  Bool minimizeCGOn; // Flag TRUE-> CG minimization active
954  Bool minVerbose; // Flag TRUE-> print extra minimization data
955  BigReal minTinyStep; // Minimization parameter
956  BigReal minBabyStep; // Minimization parameter
957  BigReal minLineGoal; // Minimization parameter
958  Bool minimizeOn; // Flag TRUE-> minimization active
959  BigReal maximumMove; // Maximum movement per timestep
960  // during minimization
961 
962  Bool sphericalBCOn; // Flag TRUE-> spherical boundary
963  // conditions are active
964  zVector sphericalCenter; // Center specified by user
965  BigReal sphericalBCk1; // First force constant for
966  // spherical BC
967  BigReal sphericalBCk2; // Second force constant for
968  // spherical BC
969  BigReal sphericalBCr1; // First radius for spherical BC
970  BigReal sphericalBCr2; // Second radius for spherical BC
971  int sphericalBCexp1; // First radius for spherical BC
972  int sphericalBCexp2; // Second radius for spherical BC
973 
974  Bool cylindricalBCOn; // Flag TRUE->cylindrical boundary
975  // conditions are active
977  char cylindricalBCAxis; // 'x', 'y', or 'z'
986 
987  Bool eFieldOn; // Should a electric field be applied
988  Bool eFieldNormalized; // Is eField vector scaled by cell basis vectors
989  zVector eField; // Electric field vector to be applied
990  BigReal eFieldFreq; // Frequency of the electric field
991  BigReal eFieldPhase; // Phase phi, cos(w*t-phi*PI/180)
992 
993  Bool stirOn; // Should a stirring torque be applied
994  char stirFilename[NAMD_FILENAME_BUFFER_SIZE]; // Stirring filename (atoms marked)
995  //do the below two even needed to be defined?
996  BigReal stirStartingTheta; // Stir starting theta offset
997  BigReal stirVel; // Stir angular velocity
998  BigReal stirK; // Stir force harmonic spring constant
999  zVector stirAxis; // Direction of stir axis
1000  zVector stirPivot; // Pivot point of stir axis
1001 
1002  Bool extraBondsOn; // read extra bonded forces
1003  Bool extraBondsCosAngles; // extra angles are cosine-based
1004  Bool extraBondsCosAnglesSetByUser; // did the user set this explicitly
1005 
1006  Bool consForceOn; // Should constant force be applied
1009 
1010  int computeEnergies; // Number of timesteps between energey evaluations
1011  int outputEnergies; // Number of timesteps between energy
1012  // outputs
1013 
1014  int outputEnergiesPrecision; // Precision of energy outputs
1015 
1016  int outputMomenta; // Number of timesteps between momentum
1017  // outputs
1018 
1019  int outputTiming; // Number of timesteps between timing
1020  // outputs
1021 
1022  Bool outputPerformance; // Calculate performance statistics?
1023  // Running stats for performance
1024  // printed on each outputTiming step
1025 
1026  int benchmarkTime; // Terminate simulation after running
1027  // for specified number of seconds,
1028  // intended for benchmarking
1029 
1030  int outputCudaTiming; // Number of timesteps between timing
1031  // outputs of CUDA code
1032 
1033  int outputPressure; // Number of timesteps between pressure
1034  // tensor outputs
1035 
1036  Bool mergeCrossterms; // Merge crossterm energy w/ dihedrals
1037 
1038  int firstTimestep; // Starting timestep. Will be 0 unless
1039  // restarting a simulation
1040 
1041  MTSChoices MTSAlgorithm; // What multiple timestep algorithm
1042  // to use
1043 
1044  int longSplitting; // What electrostatic splitting
1045  // to use
1046 
1047  Bool ignoreMass; // Mass < 3.5 does not indicate hydrogen, etc.
1048 
1049  int splitPatch; // How are patches determined?
1050  BigReal hgroupCutoff; // what is the added hydrogen margin?
1051 
1052  int mollyOn; // mollify long range forces?
1053  BigReal mollyTol; // error tolerance for molly
1054  int mollyIter; // max number of iterations for molly
1055 
1056  int rigidBonds; // what type of rigid bonds to hydrogens
1057  // none, all, or only water
1058 
1059  BigReal rigidTol; // error tolerance for rigid bonds
1060  int rigidIter; // Number of NR iterations
1061  int rigidDie; // die if rigidTol not achieved
1062 
1063  Bool useSettle; // Use SETTLE; requires rigid waters
1064 
1065  Bool testOn; // Do tests rather than simulation
1066  Bool commOnly; // Don't do any force evaluations
1067  Bool statsOn; // Don't do any force evaluations
1068 
1069  int totalAtoms; // Total Number of atoms in simulation
1070  int maxSelfPart; // maximum number of self partitions
1071  // that a patch can be split into
1072  int maxPairPart; // maximum number of pair partitions
1073  // that a patch can be split into
1074  int numAtomsSelf; // maximum number of atoms in a single
1075  // self-compute
1076  int numAtomsSelf2; // maximum number of atoms in a pair compute
1077  // in the presence of twoAwayX,Y,Z options
1078  int numAtomsPair; // maximum number of atoms in a single
1079  // pair-compute
1080  int numAtomsPair2; // maximum number of atoms in a single
1081  // pair-compute
1082  int minAtomsPerPatch; // minimum average atoms per patch
1083  // (may create larger patches)
1084  int emptyPatchLoad; // atoms worth of load generated by empty patch
1085  // (added to atom count during node assignment)
1086  int maxExclusionFlags; // maximum size of exclusion check list
1087  // for any given atom
1088  Bool outputPatchDetails; // print number of atoms per patch
1089  Bool staticAtomAssignment; // never migrate atoms
1090  Bool replicaUniformPatchGrids; // same patch grid size on all replicas
1091 
1092  //
1093  // hydrogen bond simulation parameters
1094  //
1095 
1096  // should the hydrogen bond term be used? If FALSE, all other
1097  // hydrogen bond parameters are unnecessary in simulation.
1099 
1100  // should the antecedent atom be used in the calculation of hbonds?
1102 
1103  // exponents used in hydrogen bond energy function:
1104  // aaAngleExp = exp for H-A-AA angle term (n)
1105  // haAngleExp = exp for D-H-A angle term (m)
1106  // distAttExp = exp for attractive A-D distance term (j)
1107  // distRepExp = exp for repulsive A-D distance term (i)
1109 
1110  // cutoff D-H-A angle, and on/off angles for switch fcn (in degrees)
1112 
1113  // cutoff distance for D-A separation in hbonds (in Angstroms), and
1114  // on/off distances for hbond radial term switching function
1116 
1117  // IMD parameters
1118  int IMDon; // enable IMD
1119  int IMDversion; // version of IMD protocol
1120  int IMDport; // port on which to listen for connections
1121  int IMDfreq; // frequency at which coordinates will be available
1122  int IMDwait; // if true, pause the simulation when there is no
1123  // connection
1124  int IMDignore; // IMD connection does not influence simulation
1125  // only sends coordinates and energies to VMD
1126  int IMDignoreForces; // Only the Forces are ignored. Finish, Pause and Resume are enabled
1127  IMDSessionInfo IMDsendsettings = {0,0,1,1,0,0,1}; // information about the current session for sending via IMD
1128  // default value is for v2.0
1129 
1130 
1131  // AMBER options
1132  Bool amberOn; // FLAG TRUE-> amber force field is used
1133  Bool oldParmReader; // FLAG TRUE -> use the old Amber parm/parm7 reader
1134  Bool readExclusions; // FLAG TRUE-> Read exclusions from parm file
1135  BigReal vdwscale14; // Scaling factor for 1-4 VDW interactions
1136 
1137  // GROMACS options
1138  Bool gromacsOn; // FLAG TRUE -> gromacs-style force field is used
1139 
1140  // OPLS options
1141  Bool vdwGeometricSigma; // Lennard-J sigma uses geometric mean
1142 
1143  // ScriptTcl argument passing
1151  char scriptStringArg1[128];
1152  char scriptStringArg2[128];
1153 
1156 
1158 
1161 
1162  //default value is -1
1165 
1167 
1168 
1169  //fields needed for Parallel IO Input
1173  char *binVelFile;
1174  char *binRefFile;
1175 
1176  //fields needed for Parallel IO Output
1179 
1180  char computeMapFilename[NAMD_FILENAME_BUFFER_SIZE]; // store compute map
1183 
1184  // MIC-specific parameters
1192 
1193  // AVX-512 Tiles optimizations
1195 
1199 
1200 
1201 public:
1202 
1203  SimParameters() : mgridforcelist(), parseopts(0) {};
1204  SimParameters(ConfigList *c, char *&cwd) : mgridforcelist(), parseopts(0) {
1205  initialize_config_data(c,cwd);
1206  };
1208 
1209  void initialize_config_data(ConfigList *, char *&cwd);
1210  // Initialize SimParameters data
1211  // from the ConfigList object
1212  void send_SimParameters(MOStream *);
1213  // Used by the master process
1214  // to send the paramters to
1215  // the other processors
1217  // Used by the other processors
1218  // to receive the data from the
1219  // master process
1220  void scriptSet(const char *, const char *);
1221  // Set parameters at run time
1222  void close_dcdfile(); // *** implemented in Output.C ***
1223  void close_veldcdfile(); // *** implemented in Output.C ***
1224  static void nonbonded_select();
1225  static void pme_select();
1226 
1231 
1233  return (fullElectFrequency != 1) || (nonbondedFrequency != 1);
1234  }
1235 
1236  char* getfromparseopts(const char* name, char *outbuf);
1237  int istrueinparseopts(const char* name);
1238  int issetinparseopts(const char* name);
1239 
1240  void readExtendedSystem(const char *filename, Lattice *latptr=0);
1241 private:
1242  ParseOptions *parseopts;
1243 
1244  void config_parser(ParseOptions &opts);
1245 
1246  void config_parser_basic(ParseOptions &opts);
1247  void config_parser_fileio(ParseOptions &opts);
1248  void config_parser_fullelect(ParseOptions &opts);
1249  void config_parser_methods(ParseOptions &opts);
1250  void config_parser_constraints(ParseOptions &opts);
1251 #ifdef OPENATOM_VERSION
1252  void config_parser_openatom(ParseOptions &opts);
1253 #endif //OPENATOM_VERSION
1254  /* BEGIN gf */
1255  void config_parser_gridforce(ParseOptions &opts);
1256  /* END gf */
1257  void config_parser_dcd_selections(ParseOptions &opts);
1258  void config_parser_movdrag(ParseOptions &opts);
1259  void config_parser_rotdrag(ParseOptions &opts);
1260  void config_parser_constorque(ParseOptions &opts);
1261  void config_parser_boundary(ParseOptions &opts);
1262  void config_parser_misc(ParseOptions &opts);
1263  void config_parser_mgridforce(ParseOptions &opts);
1264  void config_parser_group_restraints(ParseOptions &opts);
1265  void parse_mgrid_string_param(ConfigList *config,
1266  const char *fieldname, char** dest);
1267  void parse_mgrid_params(ConfigList *config);
1268  void print_mgrid_params();
1269  void parse_group_restraints_params(ConfigList *config);
1270 
1271  void check_config(ParseOptions &opts, ConfigList *config, char *&cwd);
1272 
1273  void print_config(ParseOptions &opts, ConfigList *config, char *&cwd);
1274 
1275  void create_output_directories(const char *dirname);
1276 
1277  int fmaFrequency; // outdated parameter name
1278  char loadBalancer[64]; // Load balancer
1279  char loadStrategy[64]; // Load balancing strategy
1280 
1281 };
1282 
1283 #endif
1284 
char qmLSSResname[5]
BigReal adaptTempTmax
int adaptTempRestartFreq
zVector sphericalCenter
Bool accelMDGresetVaftercmd
char movDragVelFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal sphericalBCr2
BigReal berendsenPressureCompressibility
Bool berendsenPressureOn
BigReal scriptArg3
char symmetryFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal TMDInitialRMSD
BigReal hgroupCutoff
char scriptStringArg1[128]
BigReal berendsenPressureRelaxationTime
BigReal soluteScalingFactorCharge
BigReal cylindricalBCl1
char qmChrgModeS[16]
char qmPCSchemeS[16]
Bool simulateInitialMapping
Bool fixedAtomsForceOutput
BigReal scriptArg2
int isSendSpanningTreeUnset()
char extCoordFilename[NAMD_FILENAME_BUFFER_SIZE]
char scriptStringArg2[128]
BigReal minTinyStep
BigReal accelMDE
BigReal monteCarloAcceptanceRate
BigReal cylindricalBCl2
BigReal adaptTempCgamma
int pressureProfileEwaldX
BigReal ldbRelativeGrainsize
int istrueinparseopts(const char *name)
BigReal getBondLambda(const BigReal) const
BigReal pairlistTrigger
void receive_SimParameters(MIStream *)
BigReal solvent_dielectric
zVector rotConsAxis
char qmBondSchemeS[16]
BigReal langevinPistonTemp
BigReal dhaOffAngle
BigReal accelMDLastStep
void close_veldcdfile()
Definition: Output.C:880
BigReal adaptTempTmin
Bool monteCarloPressureOn
Bool extraBondsCosAngles
int movingAverageWindowSize
BigReal LJPMEEwaldCoefficient
Bool globalMasterScaleByFrequency
BigReal scale14alt
int isRecvSpanningTreeUnset()
Bool ignoreExclusionChecksum
char extForcesCommand[NAMD_FILENAME_BUFFER_SIZE]
zVector sphericalConstrCenter
char rotDragAxisFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal alchElecLambdaStart
char adaptTempInFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal eFieldFreq
zVector monteCarloMaxVolume
zVector gridforceVOffset
char stirFilename[NAMD_FILENAME_BUFFER_SIZE]
Bool useGPUAtomMigration
char consTorquePivotFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal scriptArg5
BigReal multigratorPressureTarget
BigReal nonbondedScaling
int accelMDGEquiPrepSteps
Bool CUDASOAintegrateMode
char * FFTWWisdomString
Bool excludeFromPressure
BigReal constraintScaling
Bool selectConstraintsOn
char velDcdFilename[NAMD_FILENAME_BUFFER_SIZE]
BigReal sphericalBCk1
Bool ldbUnloadOutputPEs
int fullDispersionFrequency
BigReal surfaceTensionTarget
BigReal COLDTemp
BigReal reassignTemp
float Real
Definition: common.h:118
int symmetryLastFullStep
Bool movingConstraintsOn
BigReal gbis_gamma
char * tclBCScript
char qmCSMDFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal dhaOnAngle
BigReal alchLambda2
BigReal limitDist
BigReal alchBondLambdaEnd
BigReal accelMDTE
char consTorqueFile[NAMD_FILENAME_BUFFER_SIZE]
void scriptSet(const char *, const char *)
BigReal rotDragGlobVel
zVector cellBasisVector1
char extForceFilename[NAMD_FILENAME_BUFFER_SIZE]
BigReal minLineGoal
BigReal minBabyStep
char qmPrepProc[NAMD_FILENAME_BUFFER_SIZE]
int berendsenPressureFreq
BigReal cylindricalBCk2
BigReal cylindricalBCk1
int monteCarloAdjustmentFreq
Bool langevinPistonBarrier
char auxFilename[NAMD_FILENAME_BUFFER_SIZE]
Bool globalMasterStaleForces
zVector rotConsPivot
BigReal accelMDalpha
zVector cellBasisVector3
#define NAMD_FILENAME_BUFFER_SIZE
Definition: common.h:45
int pressureProfileSlabs
BigReal getElecLambda(const BigReal) const
BigReal FMMPadding
char computeMapFilename[NAMD_FILENAME_BUFFER_SIZE]
BigReal alchLambda
Bool pairInteractionOn
char qmExecPath[NAMD_FILENAME_BUFFER_SIZE]
char symmetryMatrixFile[NAMD_FILENAME_BUFFER_SIZE]
void close_dcdfile()
Definition: Output.C:864
char outputFilename[NAMD_FILENAME_BUFFER_SIZE]
BigReal switchingDist
BigReal coulomb_radius_offset
void initialize_config_data(ConfigList *, char *&cwd)
Bool langevinOnAtStartup
BigReal drudeNbtholeCut
Bool symmetryScaleForces
Bool sphericalConstraintsOn
BigReal vdwscale14
Bool useDeviceMigration
BigReal eFieldPhase
char tabulatedEnergiesFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal accelMDTalpha
WaterModel watmodel
char tclBCArgs[NAMD_FILENAME_BUFFER_SIZE]
BigReal alchRepLambdaEnd
BigReal langevinPistonDecay
char qmPCSwitchTypeS[16]
char qmSecProc[NAMD_FILENAME_BUFFER_SIZE]
BigReal rotConsVel
Bool useGPUNonbondedForceTable
BigReal ldbHomeBackgroundScaling
BigReal alchLambdaIDWS
zVector stirAxis
BigReal alpha_max
Bool langevin_useBAOAB
BigReal berendsenPressureTarget
char qmLSSModeS[16]
Bool groupRestraintsOn
zVector cylindricalCenter
BigReal LJPMEGridSpacing
BigReal stochRescalePeriod
uint32_t uint32
Definition: common.h:42
char cylindricalBCAxis
zVector gridforceScale
int accelMDGcMDPrepSteps
BigReal rescaleTemp
BigReal surface_tension
BigReal ion_concentration
BigReal PMEGridSpacing
bool isMultiTimeStepping()
BigReal accelMDFirstStep
Bool adaptTempLangevin
BigReal soluteScalingFactor
char TMDFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal dhaCutoffAngle
char rotDragPivotFile[NAMD_FILENAME_BUFFER_SIZE]
zVector strainRate2
BigReal langevinDamping
GoChoices goMethod
char dcdFilename[NAMD_FILENAME_BUFFER_SIZE]
Bool gridforcechecksize
Bool staticAtomAssignment
int Bool
Definition: common.h:142
Bool replicaUniformPatchGrids
BigReal drudeBondLen
static void pme_select()
char adaptTempRestartFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal gbis_beta
BigReal langevinTemp
MTSChoices MTSAlgorithm
BigReal scriptArg1
int monteCarloPressureFreq
BigReal loweAndersenRate
void readExtendedSystem(const char *filename, Lattice *latptr=0)
BigReal alpha_cutoff
BigReal PMEEwaldCoefficient
Bool useCUDANonbondedForceTable
BigReal fmaTheta
char alchOutFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal MSMGridSpacing
int multigratorNoseHooverChainLength
BigReal multigratorTemperatureTarget
BigReal multigratorPressureRelaxationTime
int ExclusionSettings
Definition: SimParameters.h:29
BigReal tableMaxDist
BigReal reassignIncr
char movDragFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal cylindricalBCr2
BigReal soluteScalingFactorVdw
zVector strainRate
BigReal sphericalBCr1
Bool pressureProfileEwaldOn
Bool LJPMESerialRealSpaceOn
BigReal alchVdwLambdaEnd
BigReal daCutoffDist
int GoChoices
Definition: SimParameters.h:87
BigReal multigratorTemperatureRelaxationTime
GroupRestraintList groupRestraints
BigReal ldbBackgroundScaling
Bool vdwForceSwitching
BigReal adaptTempDt
MGridforceParamsList mgridforcelist
BigReal drudeBondConst
static void nonbonded_select()
BigReal consTorqueGlobVal
BigReal accelMDGSigma0P
IMDSessionInfo IMDsendsettings
char restartFilename[NAMD_FILENAME_BUFFER_SIZE]
char symmetrykfile[NAMD_FILENAME_BUFFER_SIZE]
BigReal langevinPistonPeriod
BigReal stirStartingTheta
int pressureProfileEwaldY
Bool qmMOPACAddConfigChrg
int isRecvSpanningTreeOn()
char * getfromparseopts(const char *name, char *outbuf)
char forceDcdFilename[NAMD_FILENAME_BUFFER_SIZE]
char qmSoftware[128]
SimParameters(ConfigList *c, char *&cwd)
char consTorqueValFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal monteCarloTemp
BigReal loweAndersenCutoff
Bool monteCarloPressureOnAtStartup
Bool alchFepOnAtStartup
unsigned int randomSeed
BigReal alchTemp
zVector cellBasisVector2
BigReal initialTemp
int pressureProfileAtomTypes
BigReal scriptArg4
char accelMDGRestartFile[NAMD_FILENAME_BUFFER_SIZE]
int MTSChoices
Definition: SimParameters.h:50
BigReal getCurrentLambda2(const int) const
Bool GPUresidentSingleProcessMode
Bool extraBondsCosAnglesSetByUser
BigReal getCurrentLambda(const int) const
Bool qmBondColumnDefined
BigReal getLambdaDelta(void) const
char consTorqueAxisFile[NAMD_FILENAME_BUFFER_SIZE]
Bool langevinGammasDiffer
BigReal TMDFinalRMSD
BigReal movDragGlobVel
BigReal symmetryk
BigReal ldbPMEBackgroundScaling
Bool alchThermIntOnAtStartup
Bool pressureProfileOn
int symmetryFirstFullStep
char rotDragVelFile[NAMD_FILENAME_BUFFER_SIZE]
char SMDFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal PMETolerance
BigReal tCoupleTemp
BigReal LJPMETolerance
BigReal pairlistDist
Bool pairInteractionSelf
char rotDragFile[NAMD_FILENAME_BUFFER_SIZE]
int multigratorPressureFreq
BigReal getVdwLambda(const BigReal) const
BigReal patchDimension
BigReal getRepLambda(const BigReal) const
BigReal maximumMove
Bool qmParamPDBDefined
BigReal COLDRate
BigReal gbis_delta
int pairInteractionGroup2
BigReal pairlistShrink
BigReal dielectric
char qmBondValueTypeS[16]
size_t alchGetNumOfPMEGrids() const
BigReal adaptTempAutoDt
int isSendSpanningTreeOn()
BigReal cylindricalBCr1
char xstFilename[NAMD_FILENAME_BUFFER_SIZE]
BigReal pairlistGrow
char FFTWWisdomFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal accelMDGSigma0D
BigReal reassignHold
Bool tabulatedEnergies
BigReal MSMPadding
char TMDFile2[NAMD_FILENAME_BUFFER_SIZE]
BigReal monteCarloPressureTarget
int pairInteractionGroup1
ExclusionSettings exclude
char qmParamPDB[NAMD_FILENAME_BUFFER_SIZE]
BigReal drudeDamping
char crashFilename[NAMD_FILENAME_BUFFER_SIZE]
char qmBaseDir[NAMD_FILENAME_BUFFER_SIZE]
char goCoordinates[NAMD_FILENAME_BUFFER_SIZE]
WaterModel
Definition: common.h:221
BigReal langevinPistonTarget
int cudaGlobalProfilingFreq
int issetinparseopts(const char *name)
char qmBondColumn[16]
BigReal sphericalBCk2
Bool noPatchesOnOutputPEs
int outputEnergiesPrecision
char tableInterpType[128]
int multigratorTemperatureFreq
zVector cellOrigin
BigReal stochRescaleTemp
int globalMasterFrequency
char consForceFile[NAMD_FILENAME_BUFFER_SIZE]
char goParameters[NAMD_FILENAME_BUFFER_SIZE]
double BigReal
Definition: common.h:123
int groupRestraintsCount
zVector movingConsVel
BigReal consForceScaling
void send_SimParameters(MOStream *)
char qmColumn[16]
BigReal loweAndersenTemp
BigReal alchVdwShiftCoeff
int pressureProfileEwaldZ
BigReal drudeTemp