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, 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
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
50 #define EGP_EXCL (1 << 0)
51 #define EGP_TABLE (1 << 1)
61 class KeyValueTreeObject;
68 //! Number of T-Coupl groups
70 //! Number of of Nose-Hoover chains per group
72 //! Number of Accelerate groups
74 //! Number of Freeze groups
76 //! Number of Energy groups
78 //! Number of degrees of freedom in a group
80 //! Coupling temperature per group
82 //! No/simple/periodic simulated annealing for each group
84 //! Number of annealing time points per group
86 //! For each group: Time points
88 //! For each group: Temperature at these times. Final temp after all intervals is ref_t
92 //! Acceleration per group
94 //! Whether the group will be frozen in each direction
96 //! Exclusions/tables of energy group pairs
100 //! Number of QM groups
106 //! Simulated temperature scaling; linear or exponential
108 //! The low temperature for simulated tempering
110 //! The high temperature for simulated tempering
112 //! The range of temperatures used for simulated tempering
118 //! The frequency for calculating dhdl
120 //! 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)
122 //! The initial number of the state
124 //! Change of lambda per time step (fraction of (0.1)
126 //! Print no, total or potential energies in dhdl
127 int edHdLPrintEnergy;
128 //! The number of foreign lambda points
130 //! The array of all lambda values
132 //! The number of neighboring lambda states to calculate the energy for in up and down directions (-1 for all)
133 int lambda_neighbors;
134 //! The first lambda to calculate energies for
136 //! The last lambda +1 to calculate energies for
138 //! Free energy soft-core parameter
140 //! Lambda power for soft-core interactions
142 //! R power for soft-core interactions
144 //! Free energy soft-core sigma when c6 or c12=0
146 //! Free energy soft-core sigma for ?????
148 //! Use softcore for the coulomb portion as well (default FALSE)
150 //! Whether to print the dvdl term associated with this term; if it is not specified as separate, it is lumped with the FEP term
151 gmx_bool separate_dvdl[efptNR];
152 //! Whether to write a separate dhdl.xvg file note: NOT a gmx_bool, but an enum
153 int separate_dhdl_file;
154 //! Whether to calculate+write dhdl derivatives note: NOT a gmx_bool, but an enum
155 int dhdl_derivatives;
156 //! The maximum table size for the dH histogram
158 //! The spacing for the dH histogram
159 double dh_hist_spacing;
164 //! The frequency of expanded ensemble state changes
166 //! Which type of move updating do we use for lambda monte carlo (or no for none)
168 //! What move set will be we using for state space moves
170 //! The method we use to decide of we have equilibrated the weights
172 //! The minumum number of samples at each lambda for deciding whether we have reached a minimum
174 //! Wang-Landau delta at which we stop equilibrating weights
176 //! Use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating
178 //! After equil_steps steps we stop equilibrating the weights
180 //! After equil_samples total samples (steps/nstfep), we stop equilibrating the weights
182 //! Random number seed for lambda mc switches
184 //! Whether to use minumum variance weighting
186 //! The number of samples needed before kicking into minvar routine
188 //! The offset for the variance in MinVar
190 //! Range of cvalues used for BAR
192 //! Whether to print symmetrized matrices
193 gmx_bool bSymmetrizedTMatrix;
194 //! How frequently to print the transition matrices
196 //! Number of repetitions in the MC lambda jumps MRS -- VERIFY THIS
198 //! Minimum number of samples for each state before free sampling MRS -- VERIFY THIS!
199 int lmc_forced_nstart;
200 //! Distance in lambda space for the gibbs interval
202 //! Scaling factor for Wang-Landau
204 //! Ratio between largest and smallest number for freezing the weights
206 //! Starting delta for Wang-Landau
208 //! Use one over t convergence for Wang-Landau when the delta get sufficiently small
209 gmx_bool bWLoneovert;
210 //! Did we initialize the weights? TODO: REMOVE FOR 5.0, no longer needed with new logic
211 gmx_bool bInit_weights;
212 //! To override the main temperature, or define it if it's not defined
214 //! User-specified initial weights to start with
215 real* init_lambda_weights;
220 //! Rotation type for this group
222 //! Use mass-weighed positions?
224 //! Number of atoms in the group
226 //! The global atoms numbers
228 //! The reference positions
230 //! The normalized rotation vector
232 //! Rate of rotation (degree/ps)
234 //! Force constant (kJ/(mol nm^2)
236 //! Pivot point of rotation axis (nm)
238 //! Type of fit to determine actual group angle
240 //! Number of angles around the reference angle for which the rotation potential is also evaluated (for fit type 'potential' only)
242 //! Distance between two angles in degrees (for fit type 'potential' only)
244 //! Slab distance (nm)
246 //! Minimum value the gaussian must have so that the force is actually evaluated
248 //! Additive constant for radial motion2 and flexible2 potentials (nm^2)
254 //! Number of rotation groups
256 //! Output frequency for main rotation outfile
258 //! Output frequency for per-slab data
266 //! Number of interactive atoms
268 //! The global indices of the interactive atoms
274 //! Name of the swap group, e.g. NA, CL, SOL
276 //! Number of atoms in this group
278 //! The global ion group atoms numbers
280 //! Requested number of molecules of this type per compartment
281 int nmolReq[eCompNR];
286 //! Period between when a swap is attempted
288 //! Use mass-weighted positions in split group
289 gmx_bool massw_split[2];
290 /*! \brief Split cylinders defined by radius, upper and lower
291 * extension. The split cylinders define the channels and are
292 * each anchored in the center of the split group */
298 //! Coupling constant (number of swap attempt steps)
300 //! Ion counts may deviate from the requested values by +-threshold before a swap is done
302 //! Offset of the swap layer (='bulk') with respect to the compartment-defining layers
303 real bulkOffset[eCompNR];
304 //! Number of groups to be controlled
306 //! All swap groups, including split and solvent
310 struct t_inputrec // NOLINT (clang-analyzer-optin.performance.Padding)
313 explicit t_inputrec(const t_inputrec&) = delete;
314 t_inputrec& operator=(const t_inputrec&) = delete;
317 //! Integration method
319 //! Number of steps to be taken
321 //! Used in checkpointing to separate chunks
323 //! Start at a stepcount >0 (used w. convert-tpr)
325 //! Frequency of energy calc. and T/P coupl. upd.
327 //! Group or verlet cutoffs
329 //! Number of steps before pairlist is generated
331 //! Number of steps after which center of mass motion is removed
333 //! Center of mass motion removal algorithm
335 //! Number of steps after which print to logfile
337 //! Number of steps after which X is output
339 //! Number of steps after which V is output
341 //! Number of steps after which F is output
343 //! Number of steps after which energies printed
345 //! Number of steps after which compressed trj (.xtc,.tng) is output
346 int nstxout_compressed;
347 //! Initial time (ps)
351 //! Precision of x in compressed trajectory file
352 real x_compression_precision;
353 //! Requested fourier_spacing, when nk? not set
354 real fourier_spacing;
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
367 //! Normal/3D ewald, or pseudo-2D LR corrections
369 //! Epsilon for PME dipole correction
370 real epsilon_surface;
371 //! Type of combination rule in LJ-PME
372 int ljpme_combination_rule;
373 //! Type of periodic boundary conditions
375 //! Periodic molecules
377 //! Continuation run: starting state is correct (ie. constrained)
378 gmx_bool bContinuation;
379 //! Temperature coupling
381 //! Interval in steps for temperature coupling
383 //! Whether to print nose-hoover chains
384 gmx_bool bPrintNHChains;
385 //! Pressure coupling
387 //! Pressure coupling type
389 //! Interval in steps for pressure coupling
391 //! Pressure coupling time (ps)
393 //! Reference pressure (kJ/(mol nm^3))
395 //! Compressibility ((mol nm^3)/kJ)
397 //! How to scale absolute reference coordinates
398 int refcoord_scaling;
399 //! The COM of the posres atoms
401 //! The B-state COM of the posres atoms
403 //! Random seed for Andersen thermostat (obsolete)
405 //! Per atom pair energy drift tolerance (kJ/mol/ps/atom) for list buffer
407 //! Short range pairlist cut-off (nm)
409 //! Radius for test particle insertion
411 //! Type of electrostatics treatment
413 //! Modify the Coulomb interaction
414 int coulomb_modifier;
415 //! Coulomb switch range start (nm)
416 real rcoulomb_switch;
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;
425 //! Type of Van der Waals treatment
427 //! Modify the Van der Waals interaction
429 //! Van der Waals switch range start (nm)
431 //! Van der Waals cutoff (nm)
433 //! Perform Long range dispersion corrections
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
441 //! Data for the FEP state
443 //! Whether to do simulated tempering
445 //! Variables for simulated tempering
446 t_simtemp* simtempvals;
447 //! Whether expanded ensembles are used
449 //! Expanded ensemble parameters
450 t_expanded* expandedvals;
451 //! Type of distance restraining
453 //! Force constant for time averaged distance restraints
455 //! Type of weighting of pairs in one restraints
457 //! Use combination of time averaged and instantaneous violations
458 gmx_bool bDisreMixed;
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
473 //! Number of iterations for convergence of steepest descent in relax_shells
475 //! Stepsize for directional minimization in relax_shells
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
483 //! Order of the LINCS Projection Algorithm
485 //! Warn if any bond rotates more than this many degrees
487 //! Number of iterations in the final LINCS step
489 //! Use successive overrelaxation for shake
491 //! Friction coefficient for BD (amu/ps)
493 //! Random seed for SD and BD
495 //! The number of walls
497 //! The type of walls
499 //! The potentail is linear for r<=wall_r_linpot
501 //! The atom type for walls
502 int wall_atomtype[2];
503 //! Number density for walls
504 real wall_density[2];
505 //! Scaling factor for the box for Ewald
506 real wall_ewald_zfac;
508 /* COM pulling data */
509 //! Do we do COM pulling?
511 //! The data for center of mass pulling
515 //! Whether to use AWH biasing for PMF calculations
517 //! AWH biasing parameters
518 gmx::AwhParams* awhParams;
520 /* Enforced rotation data */
521 //! Whether to calculate enforced rotation potential(s)
523 //! The data for enforced rotation potentials
526 //! Whether to do ion/water position exchanges (CompEL)
528 //! Swap data structure.
531 //! Whether the tpr makes an interactive MD session possible.
533 //! Interactive molecular dynamics
536 //! Acceleration for viscosity calculation
538 //! Triclinic deformation velocities (nm/ps)
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
559 //! Whether twin-range scheme is active - always false if a valid .tpr was read
560 gmx_bool useTwinRange;
562 //! KVT object that contains input parameters converted to the new style.
563 gmx::KeyValueTreeObject* params;
565 //! KVT for storing simulation parameters that are not part of the mdp file.
566 std::unique_ptr<gmx::KeyValueTreeObject> internalParameters;
569 int ir_optimal_nstcalcenergy(const t_inputrec* ir);
571 int tcouple_min_integration_steps(int etc);
573 int ir_optimal_nsttcouple(const t_inputrec* ir);
575 int pcouple_min_integration_steps(int epc);
577 int ir_optimal_nstpcouple(const t_inputrec* ir);
579 /* Returns if the Coulomb force or potential is switched to zero */
580 gmx_bool ir_coulomb_switched(const t_inputrec* ir);
582 /* Returns if the Coulomb interactions are zero beyond the rcoulomb.
583 * Note: always returns TRUE for the Verlet cut-off scheme.
585 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir);
587 /* As ir_coulomb_is_zero_at_cutoff, but also returns TRUE for user tabulated
588 * interactions, since these might be zero beyond rcoulomb.
590 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir);
592 /* Returns if the Van der Waals force or potential is switched to zero */
593 gmx_bool ir_vdw_switched(const t_inputrec* ir);
595 /* Returns if the Van der Waals interactions are zero beyond the rvdw.
596 * Note: always returns TRUE for the Verlet cut-off scheme.
598 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir);
600 /* As ir_vdw_is_zero_at_cutoff, but also returns TRUE for user tabulated
601 * interactions, since these might be zero beyond rvdw.
603 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir);
605 /*! \brief Free memory from input record.
607 * All arrays and pointers will be freed.
609 * \param[in] ir The data structure
611 void done_inputrec(t_inputrec* ir);
613 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, gmx_bool bMDPformat);
615 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol);
617 void comp_pull_AB(FILE* fp, pull_params_t* pull, real ftol, real abstol);
620 gmx_bool inputrecDeform(const t_inputrec* ir);
622 gmx_bool inputrecDynamicBox(const t_inputrec* ir);
624 gmx_bool inputrecPreserveShape(const t_inputrec* ir);
626 gmx_bool inputrecNeedMutot(const t_inputrec* ir);
628 gmx_bool inputrecTwinRange(const t_inputrec* ir);
630 gmx_bool inputrecExclForces(const t_inputrec* ir);
632 gmx_bool inputrecNptTrotter(const t_inputrec* ir);
634 gmx_bool inputrecNvtTrotter(const t_inputrec* ir);
636 gmx_bool inputrecNphTrotter(const t_inputrec* ir);
638 /*! \brief Return true if the simulation is 2D periodic with two walls. */
639 bool inputrecPbcXY2Walls(const t_inputrec* ir);
641 /*! \brief Returns true for MD integator with T and/or P-coupling that supports
642 * calculating a conserved energy quantity.
644 * Note that dynamical integrators without T and P coupling (ie NVE)
645 * return false, i.e. the return value refers to whether there
646 * is a conserved quantity different than the total energy.
648 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir);
650 /*! \brief Returns true when temperature is coupled or constant. */
651 bool integratorHasReferenceTemperature(const t_inputrec* ir);
653 /*! \brief Return the number of bounded dimensions
655 * \param[in] ir The input record with MD parameters
656 * \return the number of dimensions in which
657 * the coordinates of the particles are bounded, starting at X.
659 int inputrec2nboundeddim(const t_inputrec* ir);
661 /*! \brief Returns the number of degrees of freedom in center of mass motion
663 * \param[in] ir The inputrec structure
664 * \return the number of degrees of freedom of the center of mass
666 int ndof_com(const t_inputrec* ir);
668 /*! \brief Returns the maximum reference temperature over all coupled groups
670 * Returns 0 for energy minimization and normal mode computation.
671 * Returns -1 for MD without temperature coupling.
673 * \param[in] ir The inputrec structure
675 real maxReferenceTemperature(const t_inputrec& ir);
677 /*! \brief Returns whether there is an Ewald surface contribution
679 bool haveEwaldSurfaceContribution(const t_inputrec& ir);
681 /*! \brief Check if calculation of the specific FEP type was requested.
683 * \param[in] ir Input record.
684 * \param[in] fepType Free-energy perturbation type to check for.
686 * \returns If the \p fepType is perturbed in this run.
688 bool haveFreeEnergyType(const t_inputrec& ir, int fepType);
690 #endif /* GMX_MDTYPES_INPUTREC_H */