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 
21 #include "MGridforceParams.h"
22 #include "GroupRestraintsParam.h"
23 
24 class ParseOptions;
25 class Communicate;
26 class ConfigList;
27 class MIStream;
28 class MOStream;
29 
30 // The class SimParameters is really just a glorified structure used to
31 // maintain the global simulation parameters. The only functions
32 // associated with the class are used to get the parameters from the
33 // ConfigList object, to send that Parameters from the master node
34 // to the other nodes, and to receive the Parameters on the other nodes.
35 
36 
37 // The following definitions are used to distinguish between possible
38 // bonded exclusion settings
39 typedef int ExclusionSettings;
40 
41 #define NONE 0
42 #define ONETWO 1
43 #define ONETHREE 2
44 #define ONEFOUR 3
45 #define SCALED14 4
46 
47 // The following definitions are used to distinguish between multiple
48 // timestep integration schemes
49 typedef int MTSChoices;
50 
51 #define NAIVE 0
52 #define VERLETI 1
53 
54 // The following definitions are used to distinuish between multiple
55 // long-short range force splittings
56 #define SHARP 0
57 #define XPLOR 1
58 #define C1 2
59 #define C2 3
60 
61 // The following definitions are used to distinguish among load
62 // balancers and their strategies
63 #define LDBAL_NONE 0
64 #define LDBAL_CENTRALIZED 1 // default
65 #define LDBAL_HYBRID 2
66 
67 #define LDBSTRAT_DEFAULT 10 // default
68 #define LDBSTRAT_COMPREHENSIVE 11
69 #define LDBSTRAT_REFINEONLY 12
70 #define LDBSTRAT_OLD 13
71 
72 // The following definitions are used to distinguish between patch-splitting
73 // strategies
74 #define SPLIT_PATCH_POSITION 0 // atom position determines patch
75 #define SPLIT_PATCH_HYDROGEN 1 // hydrogen groups are not broken up
76 
77 // The following definitions are used to distinguish the range of rigid
78 // bond calculations: none, all bonds to hydrogen, or only water
79 #define RIGID_NONE 0
80 #define RIGID_ALL 1
81 #define RIGID_WATER 2
82 
83 // Added by JLai -- The following definitions are used to distinguish
84 // the different GoMethodologies available to the Go program
85 // -- 6.3.11
86 typedef int GoChoices;
87 #define GO_MATRIX 1
88 #define GO_SPARSE 2
89 #define GO_LOWMEM 3
90 
91 // Used for controlling PME parallelization with ckloop
92 // The higher level will include all parallelization for lower ones
93 // E.g. If setting useCkLoop to 3, then xpencil's kspace, all
94 // backward ffts and send_untrans/ungrid routines will be parallelized
95 #define CKLOOP_CTRL_PME_UNGRIDCALC 6
96 #define CKLOOP_CTRL_PME_FORWARDFFT 5
97 #define CKLOOP_CTRL_PME_SENDTRANS 4
98 #define CKLOOP_CTRL_PME_KSPACE 3
99 #define CKLOOP_CTRL_PME_BACKWARDFFT 2
100 #define CKLOOP_CTRL_PME_SENDUNTRANS 1
101 
103 {
104 private:
105 public:
106 
107 // MAKE SURE THAT THIS CLASS CAN BE BIT COPIED OR YOU WILL HAVE TO
108 // ADD SPECIAL CODE TO send_SimParameters() and receive_SimParameters()
109 
110  int mshakeOn;
111  int lincsOn;
112 
113 #if defined(NAMD_NVTX_ENABLED) || defined(NAMD_CMK_TRACE_ENABLED) || defined(NAMD_ROCTX_ENABLED)
114  int beginEventPatchID;
115  int endEventPatchID;
116  int beginEventStep;
117  int endEventStep;
118 #endif
119 
120 #ifdef TIMER_COLLECTION
121  double timerBinWidth; // default 1
122 #endif
123 
124  Bool SOAintegrateOn; // use SOA integration routine for higher performance
125 
134  // this is set true to indicate dynamics is being performed
135 
137  // window size for moving averages of reductions for GPU-resident mode
138 
139  Bool nsPerDayOn; // prints ns/day instead of days/ns
140 
141 
142  Bool lonepairs; // enable lone pairs
143  WaterModel watmodel; // integer code for the water model in use
144  // choices are defined in common.h
145  Bool LJcorrection; // flag for whether water tail corrections should be used
146  Bool LJcorrectionAlt; // flag for whether alternative tail corrections should be used
147  BigReal dt; // Timestep size
148  int N; // Number of steps to be performed
149  int stepsPerCycle; // Number of timesteps per cycle
150 
151  zVector cellBasisVector1; // Basis vector for periodic cell
152  zVector cellBasisVector2; // Basis vector for periodic cell
153  zVector cellBasisVector3; // Basis vector for periodic cell
154  zVector cellOrigin; // Fixed center of periodic cell
155  Lattice lattice; // All data for periodic cell
156 
157  int nonbondedFrequency; // Number of timesteps between
158  // nonbonded evaluation
159  int fullElectFrequency; // Number of timesteps between
160  // full electrostatic evaluation
161  int fullDispersionFrequency; // Number of timesteps between
162  // full LJ dispersion evaluation
163  BigReal fmaTheta; // DPMTA theta value
164  int ldBalancer; // None, Centralized or Hybrid
165  int ldbStrategy; // What load balancing strategy to use
166  int ldbPeriod; // How often to do load balancing
167  int firstLdbStep; // What step to do the first
168  // load-balance on.
169  int lastLdbStep; // What step to do the last
170  // load-balance on.
171  int hybridGroupSize; // hybrid group size
172  BigReal ldbBackgroundScaling; // scaling factor for background load
173  BigReal ldbPMEBackgroundScaling;// scaling factor for PME background
174  BigReal ldbHomeBackgroundScaling;// scaling factor for home background
175  BigReal ldbRelativeGrainsize; // fraction of average load per compute
176 
177  int traceStartStep; //the timestep when trace is turned on, default to 3*firstLdbStep;
178  int numTraceSteps; //the number of timesteps that are traced, default to 2*ldbPeriod;
179 
180 #ifdef MEASURE_NAMD_WITH_PAPI
181  Bool papiMeasure; //default to false
182  int papiMeasureStartStep; //the timestep when to measure using PAPI, default to 3*firstLdbStep;
183  int numPapiMeasureSteps; //the number of timesteps when performance are measured with PAPI, default to 40;
184 #endif
185 
186  Bool outputMaps; //control whether to dump compute/patch map before load balancing
187  Bool simulateInitialMapping; //if true, the initial mapping during startup is dumped and exit
190  Bool disableTopology; // ignore torus information during patch placement
191  Bool verboseTopology; // print torus information during patch placement
192 
193  Bool benchTimestep; //only cares about benchmarking the timestep, so no file output to save SUs for large-scale benchmarking
194 
195  //whether to use CkLoop library to parallelize a loop in a function like OpenMP.
196  //It has multiple control levels. The higher the value is (must be positive), the more parallelization will be performed
197  //Currently, it is mainly used for PME computation. The default value is 0, meaning it is disabled
198  //Refer to macros CKLOOP_CTRL_* in this file for the ordering of different levels
199  int useCkLoop;
200 
201  int twoAwayX; // half-size patches in X dimension
202  int twoAwayY; // half-size patches in Y dimension
203  int twoAwayZ; // half-size patches in Z dimension
204  int maxPatches; // maximum patch count
205  Bool ldbUnloadPME; // unload processors doing PME
206  Bool ldbUnloadZero; // unload processor 0
207  Bool ldbUnloadOne; // unload processor 1
208  Bool ldbUnloadOutputPEs; // unload output processors
209  Bool noPatchesOnZero; // no patches on processor 0
210  Bool noPatchesOnOutputPEs; // no patches on output PEs
211  Bool noPatchesOnOne; // no patches on processor 1
212 
213  BigReal initialTemp; // Initial temperature for the
214  // simulation
215  Bool comMove; // Should the center of mass be
216  // able to move
217  Bool zeroMomentum; // remove momentum drift from PME
218  Bool zeroMomentumAlt; // alternate method for testing
219  Bool wrapWater; // Wrap water around on output
220  Bool wrapAll; // Wrap clusters around on output
221  Bool wrapNearest; // Wrap to closest image to origin
222  BigReal dielectric; // Dielectric constant
223  ExclusionSettings exclude; // What electrostatic exclusions should
224  // be made
225 
226  // scale14alt will override if scale14 is not set
227  BigReal scale14; // Scaling factor "1-4scaling"
228  // for 1-4 electrostatics
229  BigReal scale14alt; // Alternatively named sim parameter
230  // "oneFourScaling" so that this can
231  // be set using scripting language
232 
233  BigReal nonbondedScaling; // Scaling factor for nonbonded forces
234  int dcdFrequency; // How often (in timesteps) should
235  // a DCD trajectory file be updated
236  int dcdUnitCell; // Whether to write unit cell information in the DCD
237  int velDcdFrequency; // How often (in timesteps) should
238  // a velocity DCD file be updated
239  int forceDcdFrequency; // How often (in timesteps) should
240  // a force DCD file be updated
241  int xstFrequency; // How often (in timesteps) should
242  // a XST trajectory file be updated
243  int dcdSelectionOn; // enable dcd selection options
244  char auxFilename[NAMD_FILENAME_BUFFER_SIZE]; // auxilary output filename
245  char dcdFilename[NAMD_FILENAME_BUFFER_SIZE]; // DCD filename
246  char velDcdFilename[NAMD_FILENAME_BUFFER_SIZE]; // Velocity DCD filename
247  char forceDcdFilename[NAMD_FILENAME_BUFFER_SIZE]; // Force DCD filename
248  char xstFilename[NAMD_FILENAME_BUFFER_SIZE]; // Extended system trajectory filename
249  char outputFilename[NAMD_FILENAME_BUFFER_SIZE]; // Output file name. This name will
250  // have .coor appended to it
251  // for the coordinates and
252  // .vel appended to
253  // it for the velocities
254  char restartFilename[NAMD_FILENAME_BUFFER_SIZE]; // Base name of the restart file
255  // std::vector <DCDParams> dcdUserDefined; //parameters for user defined DCD lists
256  int restartFrequency; // How often (in timesteps) shoud the
257  // restart files be updated
258  Bool restartSave; // unique filenames for restart files
259  Bool restartSaveDcd; // unique filenames for DCD files
260  Bool binaryRestart; // should restart files be
261  // binary format rather than PDB
262  Bool binaryOutput; // should output files be
263  // binary format rather than PDB
264  BigReal cutoff; // Cutoff distance
265  BigReal margin; // Fudge factor on patch size
266  BigReal patchDimension; // Dimension of each side of a patch
267  // This is either cutoff+margin or
268  // pairlistDist+margin depending on
269  // whether or not switching is on
270  // or not
271  BigReal limitDist; // Distance below which nonbonded
272  // forces between atoms are limited
273  Bool switchingActive; // Flag TRUE->using switching function
274  // for electrostatics and vdw
275  Bool vdwForceSwitching; // Flag TRUE->using force switching
276  // function for vdw
277  BigReal switchingDist; // Distance at which switching
278  // becomes active
279  Bool martiniSwitching; // Flag TRUE->use Martini residue-based
280  // coarse-grain switching function
281  Bool martiniDielAllow; // Allow non-standard dielectric constant
282  // for use with Martini when dielectric != 15.0
283  BigReal pairlistDist; // Distance within which atom pairs
284  // should be added to pairlist
285  int pairlistMinProcs; // Minimum number of processors
286  // to enable pairlists
287  int usePairlists; // Derived from pairlistMinProcs
288 
289  int pairlistsPerCycle; // regenerate x times per cycle
290  BigReal pairlistShrink; // tol *= (1 - x) on regeneration
291  BigReal pairlistGrow; // tol *= (1 + x) on trigger
292  BigReal pairlistTrigger; // trigger is atom > (1 - x) * tol
293  int outputPairlists; // print pairlist warnings this often
294 
295  Bool constraintsOn; // Flag TRUE-> harmonic constraints
296  // active
297  int constraintExp; // Exponent for harmonic constraints
298 
299  /* BEGIN gf */
300  Bool gridforceOn; // Flag TRUE -> gridforce active
301  Bool gridforceVolts; // Flag TRUE -> gridforce using volts as units
302  zVector gridforceScale; // Gridforce scale factor
303  Bool gridforceContA1; // Flag TRUE -> grid continuous in A1 direction
304  Bool gridforceContA2; // Flag TRUE -> grid continuous in A2 direction
305  Bool gridforceContA3; // Flag TRUE -> grid continuous in A3 direction
306  zVector gridforceVOffset; // Gridforce potential offsets
307  Bool gridforceLite; // Flag TRUE -> use lightweight, fast, feature-poor gridforce
308  Bool gridforcechecksize; //Flag TRUE -> check if grid is larger than PBC cell dimensions
309  /* END gf */
312 
313  /* Begin Group restraints */
314  Bool groupRestraintsOn; // Turn on the group restraints
315  GroupRestraintList groupRestraints; // To store restraints parameters for each group
316  /* End Group restraints */
317 
318  //****** BEGIN selective restraints (X,Y,Z) changes
319  Bool selectConstraintsOn; // Flag TRUE-> selective restraints
320  // active
322  constrZOn; // Flag TRUE-> select which Cartesian
323  // component to restrain
324  //****** END selective restraints (X,Y,Z) changes
325 
326  // spherical constraints
329 
330  BigReal constraintScaling; // Scaling factor for constraint forces
331 
332  //****** BEGIN CHARMM/XPLOR type changes
333  Bool paraTypeXplorOn; // FLAG TRUE-> parametrs are XPLOR format (default)
334  Bool paraTypeCharmmOn; // FLAG TRUE-> parametrs are CHARMM format
335  //****** END CHARMM/XPLOR type changes
336 
337  // Ported by JLai -- JE - Go
338  Bool goGroPair; // FLAG FALSE->Explicit Gromacs pairs will be calculated
339  Bool goForcesOn; // FLAG TRUE-> Go forces will be calculated
340  char goParameters[NAMD_FILENAME_BUFFER_SIZE]; // File for Go parameters
341  char goCoordinates[NAMD_FILENAME_BUFFER_SIZE]; // File for Go structure and atom chain types
342  //JLai 6.3.11
343  GoChoices goMethod; // Integer for Go method -- 1) Matrix-Go, 3) Low-mem-Go
344  // End of port -- JL
345 
346  //****** BEGIN moving constraints changes
347  Bool movingConstraintsOn; // Flag TRUE-> moving constraints
348  // active
349  zVector movingConsVel; // Velocity of the movement, A/timestep
350  //****** END moving constraints changes
351  //****** BEGIN rotating constraints changes
352  Bool rotConstraintsOn; // Flag TRUE-> rotating constraints
353  // active
354  zVector rotConsAxis; // Axis of rotation
355  zVector rotConsPivot; // Pivot point of rotation
356  BigReal rotConsVel; // Velocity of rotation, Deg/timestep
357  //****** END rotating constraints changes
358 
359  //****** BEGIN moving drag changes
360  Bool movDragOn; // Flag TRUE-> moving drag active
361  char movDragFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file defining dragged atoms
362  // by non-zero value in the column
363  BigReal movDragGlobVel; // global drag velocity (A/step)
364  char movDragVelFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ scale moving drag
365  // velocity for each atom
366  //****** END moving drag changes
367  //****** BEGIN rotating drag changes
368  Bool rotDragOn; // Flag TRUE-> rotating drag active
369  char rotDragFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file defining dragged atoms
370  // by non-zero value in the column
371  char rotDragAxisFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ define axes for atoms;
372  char rotDragPivotFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ define pivots for atoms
373  BigReal rotDragGlobVel; // global drag velocity (deg/step)
374  char rotDragVelFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; B or O scales angular
375  // velocity for each atom
376  //****** END rotating drag changes
377  //****** BEGIN "constant" torque changes
378  Bool consTorqueOn; // Flag TRUE-> "constant" torque active
379  char consTorqueFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file defining torqued atoms
380  // by non-zero value in the column
381  char consTorqueAxisFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; XYZ define axes for atoms;
382  char consTorquePivotFile[NAMD_FILENAME_BUFFER_SIZE];// PDB file; XYZ define pivots for atoms
383  BigReal consTorqueGlobVal; // global "torque" (Kcal/(mol*A^2))
384  char consTorqueValFile[NAMD_FILENAME_BUFFER_SIZE]; // PDB file; B or O scales "torque"
385  // for each atom
386  //****** END "constant" torque changes
387 
388  //****** BEGIN SMD constraints changes
389  Bool SMDOn; // Flag TRUE-> SMD constraints active
390  BigReal SMDVel; // Velocity of the movement, A/timestep
391  zVector SMDDir; // Direction of the movement
392  BigReal SMDk; // Elastic constant for SMD
393  BigReal SMDk2; // Transverse elastic constant for SMD
394  char SMDFile[NAMD_FILENAME_BUFFER_SIZE]; // File for SMD information
395  int SMDOutputFreq; // Output frequency for SMD constr.
396  //****** END SMD constraints changes
397 
398  //****** BEGIN tabulated energy section
402  char tableInterpType[128];
405  //****** END tabulated energy section
406 
407  // TMD
414 
415  //Symmetry restraints
422 
423 
424 //fepb
425  Bool alchOnAtStartup; // Ensure that alchemy is set up properly
428  Bool alchOn; // Doing alchemical simulation?
429  Bool alchFepOn; // Doing alchemical simulation?
430  Bool singleTopology; // Using single topology setup?
431  Bool sdScaling; // Scaling S-D bond terms in single topology?
432  Bool alchThermIntOn; // Doing thermodynamic integration?
433  Bool alchWCAOn; // Using WCA decomposition for vdWs?
434  int alchMethod; // Which alchemical method to use? fep or ti
435  BigReal alchLambda; // lambda for dynamics
436  BigReal alchLambda2; // lambda for comparison
437  BigReal alchLambdaIDWS; // alternate lambda for interleaved double-wide sampling
438  int alchIDWSFreq; // freq with which lambda2 changes to lambdaIDWS
439  int alchLambdaFreq; // freq. (in steps) with which lambda changes
440  // from alchLambda to alchLambda2
441  BigReal getCurrentLambda(const int) const; // getter for changing lambda
442  BigReal getCurrentLambda2(const int) const; // getter for alternating lambda2 in IDWS
443  int setupIDWS(); // activates IDWS and sets alchIDWSFreq
444  BigReal getLambdaDelta(void) const; // getter for lambda increment
445  BigReal alchTemp; // temperature for alchemical calculation
446  int alchOutFreq; // freq. of alchemical output
447  Bool alchEnsembleAvg; //if do ensemble average for the net free energy difference
448  char alchOutFile[NAMD_FILENAME_BUFFER_SIZE]; // alchemical output filename
449  int alchEquilSteps; // # of equil. steps in the window
450  BigReal alchVdwShiftCoeff; // r2 shift coeff used for generating
451  // the alchemical altered vdW interactions
452  BigReal alchElecLambdaStart; // lambda value for starting point of
453  // electrostatic interactions of
454  // exnihilated particles. For annihilated
455  // particles the starting point is
456  // (1-alchElecLambdaStart)
457  BigReal getElecLambda(const BigReal) const; // return min[0,x/(1-elecStart)]
458  BigReal alchVdwLambdaEnd; // lambda value for endpoint of vdW
459  // interactions of exnihilated particles.
460  // For annihilated particles the endpoint is
461  // (1-alchVdwLambdaEnd)
462  BigReal getVdwLambda(const BigReal) const; // return max[1,x/vdwEnd]
463  BigReal alchRepLambdaEnd; // lambda value for endpoint of repulsive vdW
464  // interactions of exnihilated particles.
465  // For annihilated particles the endpoint is
466  // (1-alchRepLambdaEnd). This also implies the
467  // START for attractive vdW interactions.
468  BigReal getRepLambda(const BigReal) const; // return max[1,x/repEnd]
469  BigReal alchBondLambdaEnd; // lambda value for endpoint of bonded
470  // interactions involving exnihilated particles.
471  // For annihilated particles the endpoint is
472  // (1-alchBondLambdaEnd)
473  BigReal getBondLambda(const BigReal) const; // return max[1,x/bondEnd]
474  Bool alchDecouple; // alchemical decoupling rather than annihilation
475  Bool alchBondDecouple; // decouple purely alchemical bonds
476  size_t alchGetNumOfPMEGrids() const; // get the number of PME grids required by alchemy
477 //fepe
478 
479 
480  Bool lesOn; // Locally enhanced sampling?
481  int lesFactor; // local enhancement factor
482  Bool lesReduceTemp; // Reduce enhanced atom temperature?
483  Bool lesReduceMass; // Reduce enhanced atom mass?
484 
485  // REST2
499  Bool extForcesOn; // Are ext command forces present?
503 
504 
505  // Defines variables for QM/MM calculations
506  Bool qmForcesOn; // Are QM/MM command forces present?
511  char qmSoftware[128];
512  char qmChrgModeS[16];
514  char qmColumn[16];
520  int qmFormat ;
523 
525  char qmBondColumn[16];
531  char qmBondSchemeS[16] ;
534  char qmPCSwitchTypeS[16];
536  char qmPCSchemeS[16];
539 
545 
547  int qmLSSFreq ;
548  char qmLSSResname[5] ;
549  char qmLSSModeS[16];
551 
554 
556  int qmOutFreq ;
558 
559  Bool printBadContacts; //print indices of bad contacts being moved downhill
560 
561  //gbis implicit solvent parameters
562  Bool GBISOn; //do generalized born implicit solvent
564  Bool GBISserOn; //do generalized born implicit solvent serial
567  BigReal kappa; //debye screening length; k = sqrt(ion concentration mol/L ) / 0.304
569  BigReal gbis_delta; //three parameters for born radius calc
572  BigReal alpha_cutoff; //pairwise cutoff for integrating born radius
573  BigReal alpha_max; //maximum allowable born radius
574  Bool LCPOOn; //do LCPO SASA for GBSA
575  BigReal surface_tension; //surface tension (kcal/mol/Ang^2) for LCPO
576 
577  Bool drudeOn; // Perform integration of Drude oscillators?
578  Bool drudeHardWallOn; // Apply maximum Drude bond length restriction?
579  BigReal drudeTemp; // (low) temperature for freezing Drude oscillators
580  BigReal drudeDamping; // Langevin damping coefficient (1/ps)
581  // defaults to langevinDamping
582  BigReal drudeBondLen; // Length beyond which to apply quartic
583  // restraining potential to Drude bond
584  BigReal drudeBondConst; // Force constant for restraining potential
585  BigReal drudeNbtholeCut; // Radius of thole pair interaction
586 
587  Bool pairInteractionOn; // Calculate pair interactions?
588  int pairInteractionGroup1; // Interaction group 1.
589  int pairInteractionGroup2; // Interaction group 2.
590  Bool pairInteractionSelf; // Compute just within group.
591 
592  Bool cosAngles; // Can some angles be cos-based
593  Bool globalForcesOn; // Are global forces present?
594  Bool tclForcesOn; // Are Tcl forces present?
602 #ifdef NAMD_TCL
603  Bool tclIsThreaded; // Is Tcl library thread-safe?
604 #endif
605  Bool tclBCOn; // Are Tcl boundary forces present
606  char *tclBCScript; // Script defining tclBC calcforces
607  char tclBCArgs[NAMD_FILENAME_BUFFER_SIZE]; // Extra args for calcforces command
608  Bool freeEnergyOn; // Doing free energy perturbation?
609  Bool miscForcesOn; // Using misc forces?
610  Bool colvarsOn; // Using the colvars module?
611 
612  Bool fixedAtomsOn; // Are there fixed atoms?
613  Bool fixedAtomsForces; // Calculate forces anyway?
614  Bool fixedAtomsForceOutput; // Output fixed forces?
615 
616  Bool langevinOnAtStartup; // Ensure that langevin is set up properly
617  Bool langevinOn; // Flag TRUE-> langevin dynamics active
618  BigReal langevinTemp; // Temperature for Langevin dynamics
619  BigReal langevinDamping; // Damping coefficient (1/ps)
620  Bool langevinHydrogen; // Flag TRUE-> apply to hydrogens
625  Bool langevin_useBAOAB; // Flag TRUE-> use the experimental BAOAB integrator for NVT instead of the BBK one
626  // See Leimkuhler and Matthews (AMRX 2012); implemented in NAMD by CM June2012
627 
628  // BEGIN LA
629  Bool loweAndersenOn; // Flag TRUE-> Lowe-Andersen dynamics active
630  BigReal loweAndersenTemp; // Temperature for Lowe-Andersen dynamics
631  BigReal loweAndersenRate; // Collision frequency for Lowe-Andersen dynamics (1/ps)
632  BigReal loweAndersenCutoff; // Cutoff radius for Lowe-Andersen dynamics
633  // END LA
634 
635  Bool globalOn; // Flag TRUE-> use global integrator
636  Bool dihedralOn; // Flag TRUE-> dihedral dynamics active
637  Bool COLDOn; // Flag TRUE-> constrained overdamped
638  // langevin dynamics active
639  BigReal COLDRate; // Damping coefficient for COLD.
640  BigReal COLDTemp; // Temperature for COLD.
641 
642  Bool tCoupleOn; // Flag TRUE-> Temperature coupling
643  // active
644  BigReal tCoupleTemp; // Temperature for temp coupling
645 
667  int rescaleFreq; // Velocity rescale frequency
668  BigReal rescaleTemp; // Temperature to rescale to
669 
670  Bool accelMDOn; // Perform accelerated MD
671  Bool accelMDdihe; // Apply boost to the dihedral potential
672  Bool accelMDdual; // dual boost mode
673  Bool accelMDDebugOn; // Debugging accelerated MD
674  BigReal accelMDFirstStep; // First aMD step
675  BigReal accelMDLastStep; // Last aMD step
676  int accelMDOutFreq; // aMD output frequency
677  BigReal accelMDE; // aMD E
678  BigReal accelMDalpha; // aMD alpha
679  BigReal accelMDTE; // E for total potential in the dual boost mode
680  BigReal accelMDTalpha; // alpha for total potential in the dual boost mode
681 
682  Bool accelMDG; // Perform Gaussian accelMD calculation
683  int accelMDGiE; // Flag to set the mode iE in Gaussian accelMD
684  int accelMDGcMDSteps; // Number of cMD steps
685  int accelMDGEquiSteps; // Number of quilibration steps after adding boost potential
686  int accelMDGcMDPrepSteps; // Number of preparation cMD steps
687  int accelMDGEquiPrepSteps; // Number of preparation equilibration steps
688  int accelMDGStatWindow; // Number of steps to calc avg and std
689  BigReal accelMDGSigma0P; // upper limit of std of total potential
690  BigReal accelMDGSigma0D; // upper limit of std of dihedral potential
691  Bool accelMDGRestart; // Flag to set use restart file in Gaussian accelMD
692  char accelMDGRestartFile[NAMD_FILENAME_BUFFER_SIZE]; // restart file name
693  Bool accelMDGresetVaftercmd; // Flag to reset potential after first accelMDGcMDSteps steps
694 
695  /* Begin Adaptive Temperature Sampling */
696  Bool adaptTempOn; // is adaptTempOn
697  Bool adaptTempDebug; // Debuggin adaptive temperature sampling
698  int adaptTempFirstStep; // First adaptTemp step
699  int adaptTempLastStep; // Last adaptTemp step
700  int adaptTempOutFreq; // adaptTemp output frequency
701  int adaptTempFreq; // Steps between adaptTemp updates
702  BigReal adaptTempTmin; // Lower temperature bound
703  BigReal adaptTempTmax; // Upper temperature bound
704  BigReal adaptTempAutoDt; // Auto jump size. Value determines upper bound, adaotTempDt determines lower bound
705  int adaptTempBins; // Number of bins to store average energy values
706  BigReal adaptTempDt; // timestep for adaptTemp updates - only affects Temperature random walk
707  BigReal adaptTempCgamma; // Cgamma variable for adaptive bin averaging Cgamma = 0 is normal Averaging. 1 > Cgamma >= 0
708  Bool adaptTempLangevin; // Couple to Langevin Thermostat
709  Bool adaptTempRescale; // Couple to Vel. Rescaling
710  char adaptTempInFile[NAMD_FILENAME_BUFFER_SIZE]; // Restart information for adaptTemp to read
711  char adaptTempRestartFile[NAMD_FILENAME_BUFFER_SIZE]; // File to write restart information
712  int adaptTempRestartFreq; // Frequency of writing restart output
713  Bool adaptTempRandom; // Do we assign random temperatures when we step out of [Tmin,Tmax]?
714  /* End Adaptive Temperature Sampling */
715 
716  int reassignFreq; // Velocity reassignment frequency
717  BigReal reassignTemp; // Temperature to reassign to
718  BigReal reassignIncr; // Added to reassignTemp each time
719  BigReal reassignHold; // Hold reassignTemp at this value
720 
721  Bool useGroupPressure; // Use group rather than atomic
722  // quantities for pressure calc
723 
724  Bool excludeFromPressure; // Flag TRUE-> some atoms not rescaled
725 
726  Bool useFlexibleCell; // Use anisotropic cell fluctuations
727  Bool useConstantArea; // x,y dimensions fixed.
728  Bool useConstantRatio; // x,y ratio fixed.
729 
730  Bool fixCellDims; // fix the cell dimensions
734 
735  Bool berendsenPressureOn; // Berendsen pressure bath
740 
741  Bool langevinPistonOn; // Langevin piston pressure control
742  Bool langevinPistonBarrier; // Turn off to extrapolate cell
747 
748  Bool monteCarloPressureOnAtStartup; // Ensure that MC pressure is set up properly
749  Bool monteCarloPressureOn; // MonteCarlo pressure control
750  BigReal monteCarloPressureTarget; // target pressure
751  BigReal monteCarloTemp; // Temperature used in Monte Carlo acceptance critiria
752  BigReal monteCarloAcceptanceRate; // target acceptance rate
753  zVector monteCarloMaxVolume; // initial maximum volume change for each dimension
754  int monteCarloAdjustmentFreq; // frequency of adjusting maximum scaling factor
755  int monteCarloPressureFreq; // frequency of applying MC pressure
756 
757  Bool multigratorOn; // Multigrator temperature and/or pressure control
765 
767 
768  Bool pressureProfileOn; // Compute lateral pressure profile?
769  int pressureProfileSlabs; // Number of slabs
770  int pressureProfileFreq; // How often to store profile data
772  Bool pressureProfileEwaldOn; // Compute Ewald contribution?
776 
778  zVector strainRate2; // off diagonal elements (xy, xz, yz)
779 
780  unsigned int randomSeed; // Seed for random number generator
781 
782  Bool FMAOn; // Flag TRUE-> FMA active
783  int FMALevels; // Number of Levels for FMA
784  int FMAMp; // Number of multipole terms for FMA
785  Bool FMAFFTOn; // FFT on/off flag for FMA
786  int FMAFFTBlock; // FFT blocking factor for FMA
787 
788  Bool fullDirectOn; // Should direct calculations of
789  // full electrostatics be performed?
790 
791  Bool MSMOn; // enable MSM (multilevel summation method)
792  // for long-range electrostatics
793 
794  int MSMQuality; // choose MSM quality 0 (low) - 3 (high), using
795  // optimal combination of approximation and splitting
796  // defaults to "low" for fastest performance
797 
798  int MSMApprox; // choose MSM approximation
799  // defaults to "cubic" (low) for fastest performance
800 
801  int MSMSplit; // choose MSM splitting function
802  // defaults to "Taylor2" (low) for fastest performance
803 
804  int MSMLevels; // select number of MSM levels
805  // default (0) adapts number of levels to the
806  // system for fastest performance
807 
808  int MSMBlockSizeX; // controls size of parallel work decomposition
809  int MSMBlockSizeY; // controls size of parallel work decomposition
810  int MSMBlockSizeZ; // controls size of parallel work decomposition
811 
812  BigReal MSMGridSpacing; // defaults to 2.5 A, best for atomic systems
813 
814  BigReal MSMPadding; // pad grid along non-periodic boundaries
815  // defaults to 2.5 A
816  // increase if atoms are drifting beyond
817  // edge of grid, which will terminate
818  // simulation prematurely
819 
820  BigReal MSMxmin; // define extent of non-periodic boundaries
826 
827  Bool MsmSerialOn; // use serial MSM solver for testing
828 
832 
833  // LJ-PME parameters
834  Bool LJPMEOn; // Flag TRUE -> LJ-PME active
835  BigReal LJPMETolerance; // LJ Direct space tolerance
836  BigReal LJPMEEwaldCoefficient; // From tolerance and cutoff
837  int LJPMEInterpOrder; // LJ-PME Order of interpolation
838  int LJPMEGridSizeX; // No. of grid points in x dim
839  int LJPMEGridSizeY; // No. of grid points in y dim
840  int LJPMEGridSizeZ; // No. of grid points in z dim
841  BigReal LJPMEGridSpacing; // Maximum spacing between points
842  //Bool LJFFTWUseWisdom; // Read/save wisdom file for FFTW
843  //char LJFFTWWisdomFile[NAMD_FILENAME_BUFFER_SIZE]; // Wisdom file name for FFTW
844  //char *LJFFTWWisdomString;
845 
846  // PME Parameters
847  Bool PMEOn; // Flag TRUE -> PME active
848  BigReal PMETolerance; // Direct space tolerance
849  BigReal PMEEwaldCoefficient; // From tolerance and cutoff
850  int PMEInterpOrder; // Order of interpolation
851  int PMEGridSizeX; // No. of grid points in x dim
852  int PMEGridSizeY; // No. of grid points in y dim
853  int PMEGridSizeZ; // No. of grid points in z dim
854  BigReal PMEGridSpacing; // Maximum spacing between points
855  int PMEProcessors; // No. of processors to use
856  int PMEMinSlices; // Min slices per PME slab
857  int PMEMinPoints; // Min points per PME pencil
858  Bool PMEBarrier; // Use barrier before sendTrans
859  int PMEPencils; // Size of pencil grid in each dim
860  int PMEPencilsX; // Size of pencil grid in X dim
861  int PMEPencilsY; // Size of pencil grid in Y dim
862  int PMEPencilsZ; // Size of pencil grid in Z dim
863  int PMEPencilsYLayout; // Y pencil layout strategy
864  int PMEPencilsXLayout; // X pencil layout strategy
865  int PMESendOrder; // Message ordering strategy
866  Bool PMEOffload; // Offload reciprocal sum to accelerator
867 
868  Bool useDPME; // Flag TRUE -> old DPME code
869 
880  int bondedGPU;
902 
903  #ifdef OPENATOM_VERSION
904  Bool openatom; // Flag TRUE -> OpenAtom QM/MM active
905  #endif // OPENATOM_VERSION
906 
907  Bool minimizeCGOn; // Flag TRUE-> CG minimization active
908  Bool minVerbose; // Flag TRUE-> print extra minimization data
909  BigReal minTinyStep; // Minimization parameter
910  BigReal minBabyStep; // Minimization parameter
911  BigReal minLineGoal; // Minimization parameter
912  Bool minimizeOn; // Flag TRUE-> minimization active
913  BigReal maximumMove; // Maximum movement per timestep
914  // during minimization
915 
916  Bool sphericalBCOn; // Flag TRUE-> spherical boundary
917  // conditions are active
918  zVector sphericalCenter; // Center specified by user
919  BigReal sphericalBCk1; // First force constant for
920  // spherical BC
921  BigReal sphericalBCk2; // Second force constant for
922  // spherical BC
923  BigReal sphericalBCr1; // First radius for spherical BC
924  BigReal sphericalBCr2; // Second radius for spherical BC
925  int sphericalBCexp1; // First radius for spherical BC
926  int sphericalBCexp2; // Second radius for spherical BC
927 
928  Bool cylindricalBCOn; // Flag TRUE->cylindrical boundary
929  // conditions are active
931  char cylindricalBCAxis; // 'x', 'y', or 'z'
940 
941  Bool eFieldOn; // Should a electric field be applied
942  Bool eFieldNormalized; // Is eField vector scaled by cell basis vectors
943  zVector eField; // Electric field vector to be applied
944  BigReal eFieldFreq; // Frequency of the electric field
945  BigReal eFieldPhase; // Phase phi, cos(w*t-phi*PI/180)
946 
947  Bool stirOn; // Should a stirring torque be applied
948  char stirFilename[NAMD_FILENAME_BUFFER_SIZE]; // Stirring filename (atoms marked)
949  //do the below two even needed to be defined?
950  BigReal stirStartingTheta; // Stir starting theta offset
951  BigReal stirVel; // Stir angular velocity
952  BigReal stirK; // Stir force harmonic spring constant
953  zVector stirAxis; // Direction of stir axis
954  zVector stirPivot; // Pivot point of stir axis
955 
956  Bool extraBondsOn; // read extra bonded forces
957  Bool extraBondsCosAngles; // extra angles are cosine-based
958  Bool extraBondsCosAnglesSetByUser; // did the user set this explicitly
959 
960  Bool consForceOn; // Should constant force be applied
963 
964  int computeEnergies; // Number of timesteps between energey evaluations
965  int outputEnergies; // Number of timesteps between energy
966  // outputs
967 
968  int outputEnergiesPrecision; // Precision of energy outputs
969 
970  int outputMomenta; // Number of timesteps between momentum
971  // outputs
972 
973  int outputTiming; // Number of timesteps between timing
974  // outputs
975 
976  Bool outputPerformance; // Calculate performance statistics?
977  // Running stats for performance
978  // printed on each outputTiming step
979 
980  int benchmarkTime; // Terminate simulation after running
981  // for specified number of seconds,
982  // intended for benchmarking
983 
984  int outputCudaTiming; // Number of timesteps between timing
985  // outputs of CUDA code
986 
987  int outputPressure; // Number of timesteps between pressure
988  // tensor outputs
989 
990  Bool mergeCrossterms; // Merge crossterm energy w/ dihedrals
991 
992  int firstTimestep; // Starting timestep. Will be 0 unless
993  // restarting a simulation
994 
995  MTSChoices MTSAlgorithm; // What multiple timestep algorithm
996  // to use
997 
998  int longSplitting; // What electrostatic splitting
999  // to use
1000 
1001  Bool ignoreMass; // Mass < 3.5 does not indicate hydrogen, etc.
1002 
1003  int splitPatch; // How are patches determined?
1004  BigReal hgroupCutoff; // what is the added hydrogen margin?
1005 
1006  int mollyOn; // mollify long range forces?
1007  BigReal mollyTol; // error tolerance for molly
1008  int mollyIter; // max number of iterations for molly
1009 
1010  int rigidBonds; // what type of rigid bonds to hydrogens
1011  // none, all, or only water
1012 
1013  BigReal rigidTol; // error tolerance for rigid bonds
1014  int rigidIter; // Number of NR iterations
1015  int rigidDie; // die if rigidTol not achieved
1016 
1017  Bool useSettle; // Use SETTLE; requires rigid waters
1018 
1019  Bool testOn; // Do tests rather than simulation
1020  Bool commOnly; // Don't do any force evaluations
1021  Bool statsOn; // Don't do any force evaluations
1022 
1023  int totalAtoms; // Total Number of atoms in simulation
1024  int maxSelfPart; // maximum number of self partitions
1025  // that a patch can be split into
1026  int maxPairPart; // maximum number of pair partitions
1027  // that a patch can be split into
1028  int numAtomsSelf; // maximum number of atoms in a single
1029  // self-compute
1030  int numAtomsSelf2; // maximum number of atoms in a pair compute
1031  // in the presence of twoAwayX,Y,Z options
1032  int numAtomsPair; // maximum number of atoms in a single
1033  // pair-compute
1034  int numAtomsPair2; // maximum number of atoms in a single
1035  // pair-compute
1036  int minAtomsPerPatch; // minimum average atoms per patch
1037  // (may create larger patches)
1038  int emptyPatchLoad; // atoms worth of load generated by empty patch
1039  // (added to atom count during node assignment)
1040  int maxExclusionFlags; // maximum size of exclusion check list
1041  // for any given atom
1042  Bool outputPatchDetails; // print number of atoms per patch
1043  Bool staticAtomAssignment; // never migrate atoms
1044  Bool replicaUniformPatchGrids; // same patch grid size on all replicas
1045 
1046  //
1047  // hydrogen bond simulation parameters
1048  //
1049 
1050  // should the hydrogen bond term be used? If FALSE, all other
1051  // hydrogen bond parameters are unnecessary in simulation.
1053 
1054  // should the antecedent atom be used in the calculation of hbonds?
1056 
1057  // exponents used in hydrogen bond energy function:
1058  // aaAngleExp = exp for H-A-AA angle term (n)
1059  // haAngleExp = exp for D-H-A angle term (m)
1060  // distAttExp = exp for attractive A-D distance term (j)
1061  // distRepExp = exp for repulsive A-D distance term (i)
1063 
1064  // cutoff D-H-A angle, and on/off angles for switch fcn (in degrees)
1066 
1067  // cutoff distance for D-A separation in hbonds (in Angstroms), and
1068  // on/off distances for hbond radial term switching function
1070 
1071  // IMD parameters
1072  int IMDon; // enable IMD
1073  int IMDport; // port on which to listen for connections
1074  int IMDfreq; // frequency at which coordinates will be available
1075  int IMDwait; // if true, pause the simulation when there is no
1076  // connection
1077  int IMDignore; // IMD connection does not influence simulation
1078  // only sends coordinates and energies to VMD
1079  int IMDignoreForces; // Only the Forces are ignored. Finish, Pause and Resume are enabled
1080 
1081 
1082  // AMBER options
1083  Bool amberOn; // FLAG TRUE-> amber force field is used
1084  Bool oldParmReader; // FLAG TRUE -> use the old Amber parm/parm7 reader
1085  Bool readExclusions; // FLAG TRUE-> Read exclusions from parm file
1086  BigReal vdwscale14; // Scaling factor for 1-4 VDW interactions
1087 
1088  // GROMACS options
1089  Bool gromacsOn; // FLAG TRUE -> gromacs-style force field is used
1090 
1091  // OPLS options
1092  Bool vdwGeometricSigma; // Lennard-J sigma uses geometric mean
1093 
1094  // ScriptTcl argument passing
1102  char scriptStringArg1[128];
1103  char scriptStringArg2[128];
1104 
1107 
1109 
1112 
1113  //default value is -1
1116 
1118 
1119 
1120  //fields needed for Parallel IO Input
1124  char *binVelFile;
1125  char *binRefFile;
1126 
1127  //fields needed for Parallel IO Output
1130 
1131  char computeMapFilename[NAMD_FILENAME_BUFFER_SIZE]; // store compute map
1134 
1135  // MIC-specific parameters
1143 
1144  // AVX-512 Tiles optimizations
1146 
1150 
1151 
1152 public:
1153 
1154  SimParameters() : mgridforcelist(), parseopts(0) {};
1155  SimParameters(ConfigList *c, char *&cwd) : mgridforcelist(), parseopts(0) {
1156  initialize_config_data(c,cwd);
1157  };
1159 
1160  void initialize_config_data(ConfigList *, char *&cwd);
1161  // Initialize SimParameters data
1162  // from the ConfigList object
1163  void send_SimParameters(MOStream *);
1164  // Used by the master process
1165  // to send the paramters to
1166  // the other processors
1168  // Used by the other processors
1169  // to receive the data from the
1170  // master process
1171  void scriptSet(const char *, const char *);
1172  // Set parameters at run time
1173  void close_dcdfile(); // *** implemented in Output.C ***
1174  void close_veldcdfile(); // *** implemented in Output.C ***
1175  static void nonbonded_select();
1176  static void pme_select();
1177 
1182 
1184  return (fullElectFrequency != 1) || (nonbondedFrequency != 1);
1185  }
1186 
1187  char* getfromparseopts(const char* name, char *outbuf);
1188  int istrueinparseopts(const char* name);
1189  int issetinparseopts(const char* name);
1190 
1191  void readExtendedSystem(const char *filename, Lattice *latptr=0);
1192 private:
1193  ParseOptions *parseopts;
1194 
1195  void config_parser(ParseOptions &opts);
1196 
1197  void config_parser_basic(ParseOptions &opts);
1198  void config_parser_fileio(ParseOptions &opts);
1199  void config_parser_fullelect(ParseOptions &opts);
1200  void config_parser_methods(ParseOptions &opts);
1201  void config_parser_constraints(ParseOptions &opts);
1202 #ifdef OPENATOM_VERSION
1203  void config_parser_openatom(ParseOptions &opts);
1204 #endif //OPENATOM_VERSION
1205  /* BEGIN gf */
1206  void config_parser_gridforce(ParseOptions &opts);
1207  /* END gf */
1208  void config_parser_dcd_selections(ParseOptions &opts);
1209  void config_parser_movdrag(ParseOptions &opts);
1210  void config_parser_rotdrag(ParseOptions &opts);
1211  void config_parser_constorque(ParseOptions &opts);
1212  void config_parser_boundary(ParseOptions &opts);
1213  void config_parser_misc(ParseOptions &opts);
1214  void config_parser_mgridforce(ParseOptions &opts);
1215  void config_parser_group_restraints(ParseOptions &opts);
1216  void parse_mgrid_string_param(ConfigList *config,
1217  const char *fieldname, char** dest);
1218  void parse_mgrid_params(ConfigList *config);
1219  void print_mgrid_params();
1220  void parse_group_restraints_params(ConfigList *config);
1221 
1222  void check_config(ParseOptions &opts, ConfigList *config, char *&cwd);
1223 
1224  void print_config(ParseOptions &opts, ConfigList *config, char *&cwd);
1225 
1226  void create_output_directories(const char *dirname);
1227 
1228  int fmaFrequency; // outdated parameter name
1229  char loadBalancer[64]; // Load balancer
1230  char loadStrategy[64]; // Load balancing strategy
1231 
1232 };
1233 
1234 #endif
1235 
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:829
BigReal adaptTempTmin
Bool monteCarloPressureOn
Bool extraBondsCosAngles
int movingAverageWindowSize
BigReal LJPMEEwaldCoefficient
Bool globalMasterScaleByFrequency
BigReal scale14alt
int isRecvSpanningTreeUnset()
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
zVector stirPivot
Bool outputPerformance
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:813
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
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:28
BigReal tableMaxDist
BigReal reassignIncr
char movDragFile[NAMD_FILENAME_BUFFER_SIZE]
BigReal cylindricalBCr2
BigReal soluteScalingFactorVdw
zVector strainRate
BigReal sphericalBCr1
Bool pressureProfileEwaldOn
BigReal alchVdwLambdaEnd
BigReal daCutoffDist
int GoChoices
Definition: SimParameters.h:86
BigReal multigratorTemperatureRelaxationTime
GroupRestraintList groupRestraints
BigReal ldbBackgroundScaling
Bool vdwForceSwitching
BigReal adaptTempDt
MGridforceParamsList mgridforcelist
BigReal drudeBondConst
static void nonbonded_select()
BigReal consTorqueGlobVal
BigReal accelMDGSigma0P
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:49
BigReal getCurrentLambda2(const int) const
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 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
zVector movingConsVel
BigReal consForceScaling
void send_SimParameters(MOStream *)
char qmColumn[16]
BigReal loweAndersenTemp
BigReal alchVdwShiftCoeff
int pressureProfileEwaldZ
BigReal drudeTemp