2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_MDTYPES_INPUTREC_H
39 #define GMX_MDTYPES_INPUTREC_H
46 #include "gromacs/math/vectypes.h"
47 #include "gromacs/mdtypes/md_enums.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49 #include "gromacs/utility/real.h"
51 #define EGP_EXCL (1 << 0)
52 #define EGP_TABLE (1 << 1)
62 class KeyValueTreeObject;
68 //! Number of T-Coupl groups
70 //! Number of of Nose-Hoover chains per group
71 int nhchainlength = 0;
72 //! Number of Freeze groups
74 //! Number of Energy groups
76 //! Number of degrees of freedom in a group
78 //! Coupling temperature per group
79 real* ref_t = nullptr;
80 //! No/simple/periodic simulated annealing for each group
81 SimulatedAnnealing* annealing = nullptr;
82 //! Number of annealing time points per group
83 int* anneal_npoints = nullptr;
84 //! For each group: Time points
85 real** anneal_time = nullptr;
86 //! For each group: Temperature at these times. Final temp after all intervals is ref_t
87 real** anneal_temp = nullptr;
89 real* tau_t = nullptr;
90 //! Whether the group will be frozen in each direction
91 ivec* nFreeze = nullptr;
92 //! Exclusions/tables of energy group pairs
93 int* egp_flags = nullptr;
96 //! Number of QM groups
102 //! Simulated temperature scaling; linear or exponential
103 SimulatedTempering eSimTempScale = SimulatedTempering::Default;
104 //! The low temperature for simulated tempering
105 real simtemp_low = 0;
106 //! The high temperature for simulated tempering
107 real simtemp_high = 0;
108 //! The range of temperatures used for simulated tempering
109 std::vector<real> temperatures;
114 //! The frequency for calculating dhdl
116 //! Fractional value of lambda (usually will use init_fep_state, this will only be for slow growth, and for legacy free energy code. Only has a valid value if positive)
117 double init_lambda = 0;
118 //! The initial number of the state
119 int init_fep_state = 0;
120 //! Change of lambda per time step (fraction of (0.1)
121 double delta_lambda = 0;
122 //! Print no, total or potential energies in dhdl
123 FreeEnergyPrintEnergy edHdLPrintEnergy = FreeEnergyPrintEnergy::Default;
124 //! The number of foreign lambda points
126 //! The array of all lambda values
127 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, std::vector<double>> all_lambda;
128 //! The number of neighboring lambda states to calculate the energy for in up and down directions (-1 for all)
129 int lambda_neighbors = 0;
130 //! The first lambda to calculate energies for
131 int lambda_start_n = 0;
132 //! The last lambda +1 to calculate energies for
133 int lambda_stop_n = 0;
134 //! Free energy soft-core parameter
136 //! Lambda power for soft-core interactions
138 //! R power for soft-core interactions
140 //! Free energy soft-core sigma when c6 or c12=0
142 //! Free energy soft-core sigma for ?????
143 real sc_sigma_min = 0;
144 //! Use softcore for the coulomb portion as well (default FALSE)
145 bool bScCoul = false;
146 //! Whether to print the dvdl term associated with this term; if it is not specified as separate, it is lumped with the FEP term
147 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, bool> separate_dvdl;
148 //! Whether to write a separate dhdl.xvg file note: NOT a gmx_bool, but an enum
149 SeparateDhdlFile separate_dhdl_file = SeparateDhdlFile::Default;
150 //! Whether to calculate+write dhdl derivatives note: NOT a gmx_bool, but an enum
151 DhDlDerivativeCalculation dhdl_derivatives = DhDlDerivativeCalculation::Default;
152 //! The maximum table size for the dH histogram
153 int dh_hist_size = 0;
154 //! The spacing for the dH histogram
155 double dh_hist_spacing = 0;
160 //! The frequency of expanded ensemble state changes
162 //! Which type of move updating do we use for lambda monte carlo (or no for none)
163 LambdaWeightCalculation elamstats = LambdaWeightCalculation::Default;
164 //! What move set will be we using for state space moves
165 LambdaMoveCalculation elmcmove = LambdaMoveCalculation::Default;
166 //! The method we use to decide of we have equilibrated the weights
167 LambdaWeightWillReachEquilibrium elmceq = LambdaWeightWillReachEquilibrium::Default;
168 //! The minumum number of samples at each lambda for deciding whether we have reached a minimum
169 int equil_n_at_lam = 0;
170 //! Wang-Landau delta at which we stop equilibrating weights
171 real equil_wl_delta = 0;
172 //! Use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating
173 real equil_ratio = 0;
174 //! After equil_steps steps we stop equilibrating the weights
176 //! After equil_samples total samples (steps/nstfep), we stop equilibrating the weights
177 int equil_samples = 0;
178 //! Random number seed for lambda mc switches
180 //! Whether to use minumum variance weighting
182 //! The number of samples needed before kicking into minvar routine
184 //! The offset for the variance in MinVar
185 real minvar_const = 0;
186 //! Range of cvalues used for BAR
188 //! Whether to print symmetrized matrices
189 bool bSymmetrizedTMatrix = false;
190 //! How frequently to print the transition matrices
192 //! Number of repetitions in the MC lambda jumps MRS -- VERIFY THIS
194 //! Minimum number of samples for each state before free sampling MRS -- VERIFY THIS!
195 int lmc_forced_nstart = 0;
196 //! Distance in lambda space for the gibbs interval
197 int gibbsdeltalam = 0;
198 //! Scaling factor for Wang-Landau
200 //! Ratio between largest and smallest number for freezing the weights
202 //! Starting delta for Wang-Landau
203 real init_wl_delta = 0;
204 //! Use one over t convergence for Wang-Landau when the delta get sufficiently small
205 bool bWLoneovert = false;
206 //! Did we initialize the weights? TODO: REMOVE FOR 5.0, no longer needed with new logic
207 bool bInit_weights = false;
208 //! To override the main temperature, or define it if it's not defined
210 //! User-specified initial weights to start with
211 std::vector<real> init_lambda_weights;
216 //! Rotation type for this group
217 EnforcedRotationGroupType eType;
218 //! Use mass-weighed positions?
220 //! Number of atoms in the group
222 //! The global atoms numbers
224 //! The reference positions
226 //! The normalized rotation vector
228 //! Rate of rotation (degree/ps)
230 //! Force constant (kJ/(mol nm^2)
232 //! Pivot point of rotation axis (nm)
234 //! Type of fit to determine actual group angle
235 RotationGroupFitting eFittype;
236 //! Number of angles around the reference angle for which the rotation potential is also evaluated (for fit type 'potential' only)
238 //! Distance between two angles in degrees (for fit type 'potential' only)
240 //! Slab distance (nm)
242 //! Minimum value the gaussian must have so that the force is actually evaluated
244 //! Additive constant for radial motion2 and flexible2 potentials (nm^2)
250 //! Number of rotation groups
252 //! Output frequency for main rotation outfile
254 //! Output frequency for per-slab data
262 //! Number of interactive atoms
264 //! The global indices of the interactive atoms
270 //! Name of the swap group, e.g. NA, CL, SOL
272 //! Number of atoms in this group
274 //! The global ion group atoms numbers
276 //! Requested number of molecules of this type per compartment
277 gmx::EnumerationArray<Compartment, int> nmolReq;
282 //! Period between when a swap is attempted
284 //! Use mass-weighted positions in split group
286 /*! \brief Split cylinders defined by radius, upper and lower
287 * extension. The split cylinders define the channels and are
288 * each anchored in the center of the split group */
294 //! Coupling constant (number of swap attempt steps)
296 //! Ion counts may deviate from the requested values by +-threshold before a swap is done
298 //! Offset of the swap layer (='bulk') with respect to the compartment-defining layers
299 gmx::EnumerationArray<Compartment, real> bulkOffset;
300 //! Number of groups to be controlled
302 //! All swap groups, including split and solvent
306 struct t_inputrec // NOLINT (clang-analyzer-optin.performance.Padding)
309 explicit t_inputrec(const t_inputrec&) = delete;
310 t_inputrec& operator=(const t_inputrec&) = delete;
313 //! Integration method
314 IntegrationAlgorithm eI = IntegrationAlgorithm::Default;
315 //! Number of steps to be taken
317 //! Used in checkpointing to separate chunks
318 int simulation_part = 0;
319 //! Start at a stepcount >0 (used w. convert-tpr)
320 int64_t init_step = 0;
321 //! Frequency of energy calc. and T/P coupl. upd.
322 int nstcalcenergy = 0;
323 //! Group or verlet cutoffs
324 CutoffScheme cutoff_scheme = CutoffScheme::Default;
325 //! Number of steps before pairlist is generated
327 //! Number of steps after which center of mass motion is removed
329 //! Center of mass motion removal algorithm
330 ComRemovalAlgorithm comm_mode = ComRemovalAlgorithm::Default;
331 //! Number of steps after which print to logfile
333 //! Number of steps after which X is output
335 //! Number of steps after which V is output
337 //! Number of steps after which F is output
339 //! Number of steps after which energies printed
341 //! Number of steps after which compressed trj (.xtc,.tng) is output
342 int nstxout_compressed = 0;
343 //! Initial time (ps)
347 //! Whether we use multiple time stepping
349 //! The multiple time stepping levels
350 std::vector<gmx::MtsLevel> mtsLevels;
351 //! Precision of x in compressed trajectory file
352 real x_compression_precision = 0;
353 //! Requested fourier_spacing, when nk? not set
354 real fourier_spacing = 0;
355 //! Number of k vectors in x dimension for fourier methods for long range electrost.
357 //! Number of k vectors in y dimension for fourier methods for long range electrost.
359 //! Number of k vectors in z dimension for fourier methods for long range electrost.
361 //! Interpolation order for PME
363 //! Real space tolerance for Ewald, determines the real/reciprocal space relative weight
365 //! Real space tolerance for LJ-Ewald
366 real ewald_rtol_lj = 0;
367 //! Normal/3D ewald, or pseudo-2D LR corrections
368 EwaldGeometry ewald_geometry = EwaldGeometry::Default;
369 //! Epsilon for PME dipole correction
370 real epsilon_surface = 0;
371 //! Type of combination rule in LJ-PME
372 LongRangeVdW ljpme_combination_rule = LongRangeVdW::Default;
373 //! Type of periodic boundary conditions
374 PbcType pbcType = PbcType::Default;
375 //! Periodic molecules
376 bool bPeriodicMols = false;
377 //! Continuation run: starting state is correct (ie. constrained)
378 bool bContinuation = false;
379 //! Temperature coupling
380 TemperatureCoupling etc = TemperatureCoupling::Default;
381 //! Interval in steps for temperature coupling
383 //! Whether to print nose-hoover chains
384 bool bPrintNHChains = false;
385 //! Pressure coupling
386 PressureCoupling epc = PressureCoupling::Default;
387 //! Pressure coupling type
388 PressureCouplingType epct = PressureCouplingType::Default;
389 //! Interval in steps for pressure coupling
391 //! Pressure coupling time (ps)
393 //! Reference pressure (kJ/(mol nm^3))
394 tensor ref_p = { { 0 } };
395 //! Compressibility ((mol nm^3)/kJ)
396 tensor compress = { { 0 } };
397 //! How to scale absolute reference coordinates
398 RefCoordScaling refcoord_scaling = RefCoordScaling::Default;
399 //! The COM of the posres atoms
400 rvec posres_com = { 0, 0, 0 };
401 //! The B-state COM of the posres atoms
402 rvec posres_comB = { 0, 0, 0 };
403 //! Random seed for Andersen thermostat (obsolete)
404 int andersen_seed = 0;
405 //! Per atom pair energy drift tolerance (kJ/mol/ps/atom) for list buffer
406 real verletbuf_tol = 0;
407 //! Short range pairlist cut-off (nm)
409 //! Radius for test particle insertion
411 //! Type of electrostatics treatment
412 CoulombInteractionType coulombtype = CoulombInteractionType::Default;
413 //! Modify the Coulomb interaction
414 InteractionModifiers coulomb_modifier = InteractionModifiers::Default;
415 //! Coulomb switch range start (nm)
416 real rcoulomb_switch = 0;
417 //! Coulomb cutoff (nm)
419 //! Relative dielectric constant
421 //! Relative dielectric constant of the RF
423 //! Always false (no longer supported)
424 bool implicit_solvent = false;
425 //! Type of Van der Waals treatment
426 VanDerWaalsType vdwtype = VanDerWaalsType::Default;
427 //! Modify the Van der Waals interaction
428 InteractionModifiers vdw_modifier = InteractionModifiers::Default;
429 //! Van der Waals switch range start (nm)
430 real rvdw_switch = 0;
431 //! Van der Waals cutoff (nm)
433 //! Perform Long range dispersion corrections
434 DispersionCorrectionType eDispCorr = DispersionCorrectionType::Default;
435 //! Extension of the table beyond the cut-off, as well as the table length for 1-4 interac.
437 //! Tolerance for shake
439 //! Free energy calculations
440 FreeEnergyPerturbationType efep = FreeEnergyPerturbationType::Default;
441 //! Data for the FEP state
442 std::unique_ptr<t_lambda> fepvals;
443 //! Whether to do simulated tempering
444 bool bSimTemp = false;
445 //! Variables for simulated tempering
446 std::unique_ptr<t_simtemp> simtempvals;
447 //! Whether expanded ensembles are used
448 bool bExpanded = false;
449 //! Expanded ensemble parameters
450 std::unique_ptr<t_expanded> expandedvals;
451 //! Type of distance restraining
452 DistanceRestraintRefinement eDisre = DistanceRestraintRefinement::Default;
453 //! Force constant for time averaged distance restraints
455 //! Type of weighting of pairs in one restraints
456 DistanceRestraintWeighting eDisreWeighting = DistanceRestraintWeighting::Default;
457 //! Use combination of time averaged and instantaneous violations
458 bool bDisreMixed = false;
459 //! Frequency of writing pair distances to enx
461 //! Time constant for memory function in disres
463 //! Force constant for orientational restraints
465 //! Time constant for memory function in orires
467 //! Frequency of writing tr(SD) to energy output
469 //! The stepsize for updating
470 real em_stepsize = 0;
473 //! Number of iterations for convergence of steepest descent in relax_shells
475 //! Stepsize for directional minimization in relax_shells
476 real fc_stepsize = 0;
477 //! Number of steps after which a steepest descents step is done while doing cg
479 //! Number of corrections to the Hessian to keep
481 //! Type of constraint algorithm
482 ConstraintAlgorithm eConstrAlg = ConstraintAlgorithm::Default;
483 //! Order of the LINCS Projection Algorithm
485 //! Warn if any bond rotates more than this many degrees
486 real LincsWarnAngle = 0;
487 //! Number of iterations in the final LINCS step
489 //! Use successive overrelaxation for shake
490 bool bShakeSOR = false;
491 //! Friction coefficient for BD (amu/ps)
493 //! Random seed for SD and BD
495 //! The number of walls
497 //! The type of walls
498 WallType wall_type = WallType::Default;
499 //! The potentail is linear for r<=wall_r_linpot
500 real wall_r_linpot = 0;
501 //! The atom type for walls
502 int wall_atomtype[2] = { 0, 0 };
503 //! Number density for walls
504 real wall_density[2] = { 0, 0 };
505 //! Scaling factor for the box for Ewald
506 real wall_ewald_zfac = 0;
508 /* COM pulling data */
509 //! Do we do COM pulling?
511 //! The data for center of mass pulling
512 std::unique_ptr<pull_params_t> pull;
515 //! Whether to use AWH biasing for PMF calculations
517 //! AWH biasing parameters
518 std::unique_ptr<gmx::AwhParams> awhParams;
520 /* Enforced rotation data */
521 //! Whether to calculate enforced rotation potential(s)
523 //! The data for enforced rotation potentials
524 t_rot* rot = nullptr;
526 //! Whether to do ion/water position exchanges (CompEL)
527 SwapType eSwapCoords = SwapType::Default;
528 //! Swap data structure.
529 t_swapcoords* swap = nullptr;
531 //! Whether the tpr makes an interactive MD session possible.
533 //! Interactive molecular dynamics
534 t_IMD* imd = nullptr;
536 //! Acceleration for viscosity calculation
538 //! Triclinic deformation velocities (nm/ps)
539 tensor deform = { { 0 } };
540 /*! \brief User determined parameters */
553 //! QM/MM calculation
556 /* Fields for removed features go here (better caching) */
557 //! Whether AdResS is enabled - always false if a valid .tpr was read
558 bool bAdress = false;
559 //! Whether twin-range scheme is active - always false if a valid .tpr was read
560 bool useTwinRange = false;
561 //! Whether we have constant acceleration - removed in GROMACS 2022
562 bool useConstantAcceleration = false;
564 //! KVT object that contains input parameters converted to the new style.
565 gmx::KeyValueTreeObject* params = nullptr;
567 //! KVT for storing simulation parameters that are not part of the mdp file.
568 std::unique_ptr<gmx::KeyValueTreeObject> internalParameters;
571 int ir_optimal_nstcalcenergy(const t_inputrec* ir);
573 int tcouple_min_integration_steps(TemperatureCoupling etc);
575 int ir_optimal_nsttcouple(const t_inputrec* ir);
577 int pcouple_min_integration_steps(PressureCoupling epc);
579 int ir_optimal_nstpcouple(const t_inputrec* ir);
581 /* Returns if the Coulomb force or potential is switched to zero */
582 bool ir_coulomb_switched(const t_inputrec* ir);
584 /* Returns if the Coulomb interactions are zero beyond the rcoulomb.
585 * Note: always returns TRUE for the Verlet cut-off scheme.
587 bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir);
589 /* As ir_coulomb_is_zero_at_cutoff, but also returns TRUE for user tabulated
590 * interactions, since these might be zero beyond rcoulomb.
592 bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir);
594 /* Returns if the Van der Waals force or potential is switched to zero */
595 bool ir_vdw_switched(const t_inputrec* ir);
597 /* Returns if the Van der Waals interactions are zero beyond the rvdw.
598 * Note: always returns TRUE for the Verlet cut-off scheme.
600 bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir);
602 /* As ir_vdw_is_zero_at_cutoff, but also returns TRUE for user tabulated
603 * interactions, since these might be zero beyond rvdw.
605 bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir);
607 /*! \brief Free memory from input record.
609 * All arrays and pointers will be freed.
611 * \param[in] ir The data structure
613 void done_inputrec(t_inputrec* ir);
615 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, bool bMDPformat);
617 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol);
619 void comp_pull_AB(FILE* fp, const pull_params_t& pull, real ftol, real abstol);
622 bool inputrecDeform(const t_inputrec* ir);
624 bool inputrecDynamicBox(const t_inputrec* ir);
626 bool inputrecPreserveShape(const t_inputrec* ir);
628 bool inputrecNeedMutot(const t_inputrec* ir);
630 bool inputrecTwinRange(const t_inputrec* ir);
632 bool inputrecExclForces(const t_inputrec* ir);
634 bool inputrecNptTrotter(const t_inputrec* ir);
636 bool inputrecNvtTrotter(const t_inputrec* ir);
638 bool inputrecNphTrotter(const t_inputrec* ir);
640 /*! \brief Return true if the simulation is 2D periodic with two walls. */
641 bool inputrecPbcXY2Walls(const t_inputrec* ir);
643 //! \brief Return true if the simulation has frozen atoms (non-trivial freeze groups).
644 bool inputrecFrozenAtoms(const t_inputrec* ir);
646 /*! \brief Returns true for MD integator with T and/or P-coupling that supports
647 * calculating a conserved energy quantity.
649 * Note that dynamical integrators without T and P coupling (ie NVE)
650 * return false, i.e. the return value refers to whether there
651 * is a conserved quantity different than the total energy.
653 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir);
655 /*! \brief Returns true when temperature is coupled or constant. */
656 bool integratorHasReferenceTemperature(const t_inputrec* ir);
658 /*! \brief Return the number of bounded dimensions
660 * \param[in] ir The input record with MD parameters
661 * \return the number of dimensions in which
662 * the coordinates of the particles are bounded, starting at X.
664 int inputrec2nboundeddim(const t_inputrec* ir);
666 /*! \brief Returns the number of degrees of freedom in center of mass motion
668 * \param[in] ir The inputrec structure
669 * \return the number of degrees of freedom of the center of mass
671 int ndof_com(const t_inputrec* ir);
673 /*! \brief Returns the maximum reference temperature over all coupled groups
675 * Returns 0 for energy minimization and normal mode computation.
676 * Returns -1 for MD without temperature coupling.
678 * \param[in] ir The inputrec structure
680 real maxReferenceTemperature(const t_inputrec& ir);
682 /*! \brief Returns whether there is an Ewald surface contribution
684 bool haveEwaldSurfaceContribution(const t_inputrec& ir);
686 /*! \brief Check if calculation of the specific FEP type was requested.
688 * \param[in] ir Input record.
689 * \param[in] fepType Free-energy perturbation type to check for.
691 * \returns If the \p fepType is perturbed in this run.
693 bool haveFreeEnergyType(const t_inputrec& ir, int fepType);
695 #endif /* GMX_MDTYPES_INPUTREC_H */