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/basedefinitions.h"
49 #include "gromacs/utility/real.h"
51 #define EGP_EXCL (1 << 0)
52 #define EGP_TABLE (1 << 1)
62 class KeyValueTreeObject;
70 //! Number of T-Coupl groups
72 //! Number of of Nose-Hoover chains per group
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 //! Whether the group will be frozen in each direction
94 //! Exclusions/tables of energy group pairs
98 //! Number of QM groups
104 //! Simulated temperature scaling; linear or exponential
106 //! The low temperature for simulated tempering
108 //! The high temperature for simulated tempering
110 //! The range of temperatures used for simulated tempering
116 //! The frequency for calculating dhdl
118 //! 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)
120 //! The initial number of the state
122 //! Change of lambda per time step (fraction of (0.1)
124 //! Print no, total or potential energies in dhdl
125 int edHdLPrintEnergy;
126 //! The number of foreign lambda points
128 //! The array of all lambda values
130 //! The number of neighboring lambda states to calculate the energy for in up and down directions (-1 for all)
131 int lambda_neighbors;
132 //! The first lambda to calculate energies for
134 //! The last lambda +1 to calculate energies for
136 //! Free energy soft-core parameter
138 //! Lambda power for soft-core interactions
140 //! R power for soft-core interactions
142 //! Free energy soft-core sigma when c6 or c12=0
144 //! Free energy soft-core sigma for ?????
146 //! Use softcore for the coulomb portion as well (default FALSE)
148 //! Whether to print the dvdl term associated with this term; if it is not specified as separate, it is lumped with the FEP term
149 gmx_bool separate_dvdl[efptNR];
150 //! Whether to write a separate dhdl.xvg file note: NOT a gmx_bool, but an enum
151 int separate_dhdl_file;
152 //! Whether to calculate+write dhdl derivatives note: NOT a gmx_bool, but an enum
153 int dhdl_derivatives;
154 //! The maximum table size for the dH histogram
156 //! The spacing for the dH histogram
157 double dh_hist_spacing;
162 //! The frequency of expanded ensemble state changes
164 //! Which type of move updating do we use for lambda monte carlo (or no for none)
166 //! What move set will be we using for state space moves
168 //! The method we use to decide of we have equilibrated the weights
170 //! The minumum number of samples at each lambda for deciding whether we have reached a minimum
172 //! Wang-Landau delta at which we stop equilibrating weights
174 //! Use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating
176 //! After equil_steps steps we stop equilibrating the weights
178 //! After equil_samples total samples (steps/nstfep), we stop equilibrating the weights
180 //! Random number seed for lambda mc switches
182 //! Whether to use minumum variance weighting
184 //! The number of samples needed before kicking into minvar routine
186 //! The offset for the variance in MinVar
188 //! Range of cvalues used for BAR
190 //! Whether to print symmetrized matrices
191 gmx_bool bSymmetrizedTMatrix;
192 //! How frequently to print the transition matrices
194 //! Number of repetitions in the MC lambda jumps MRS -- VERIFY THIS
196 //! Minimum number of samples for each state before free sampling MRS -- VERIFY THIS!
197 int lmc_forced_nstart;
198 //! Distance in lambda space for the gibbs interval
200 //! Scaling factor for Wang-Landau
202 //! Ratio between largest and smallest number for freezing the weights
204 //! Starting delta for Wang-Landau
206 //! Use one over t convergence for Wang-Landau when the delta get sufficiently small
207 gmx_bool bWLoneovert;
208 //! Did we initialize the weights? TODO: REMOVE FOR 5.0, no longer needed with new logic
209 gmx_bool bInit_weights;
210 //! To override the main temperature, or define it if it's not defined
212 //! User-specified initial weights to start with
213 real* init_lambda_weights;
218 //! Rotation type for this group
220 //! Use mass-weighed positions?
222 //! Number of atoms in the group
224 //! The global atoms numbers
226 //! The reference positions
228 //! The normalized rotation vector
230 //! Rate of rotation (degree/ps)
232 //! Force constant (kJ/(mol nm^2)
234 //! Pivot point of rotation axis (nm)
236 //! Type of fit to determine actual group angle
238 //! Number of angles around the reference angle for which the rotation potential is also evaluated (for fit type 'potential' only)
240 //! Distance between two angles in degrees (for fit type 'potential' only)
242 //! Slab distance (nm)
244 //! Minimum value the gaussian must have so that the force is actually evaluated
246 //! Additive constant for radial motion2 and flexible2 potentials (nm^2)
252 //! Number of rotation groups
254 //! Output frequency for main rotation outfile
256 //! Output frequency for per-slab data
264 //! Number of interactive atoms
266 //! The global indices of the interactive atoms
272 //! Name of the swap group, e.g. NA, CL, SOL
274 //! Number of atoms in this group
276 //! The global ion group atoms numbers
278 //! Requested number of molecules of this type per compartment
279 int nmolReq[eCompNR];
284 //! Period between when a swap is attempted
286 //! Use mass-weighted positions in split group
287 gmx_bool massw_split[2];
288 /*! \brief Split cylinders defined by radius, upper and lower
289 * extension. The split cylinders define the channels and are
290 * each anchored in the center of the split group */
296 //! Coupling constant (number of swap attempt steps)
298 //! Ion counts may deviate from the requested values by +-threshold before a swap is done
300 //! Offset of the swap layer (='bulk') with respect to the compartment-defining layers
301 real bulkOffset[eCompNR];
302 //! Number of groups to be controlled
304 //! All swap groups, including split and solvent
308 struct t_inputrec // NOLINT (clang-analyzer-optin.performance.Padding)
311 explicit t_inputrec(const t_inputrec&) = delete;
312 t_inputrec& operator=(const t_inputrec&) = delete;
315 //! Integration method
317 //! Number of steps to be taken
319 //! Used in checkpointing to separate chunks
321 //! Start at a stepcount >0 (used w. convert-tpr)
323 //! Frequency of energy calc. and T/P coupl. upd.
325 //! Group or verlet cutoffs
327 //! Number of steps before pairlist is generated
329 //! Number of steps after which center of mass motion is removed
331 //! Center of mass motion removal algorithm
333 //! Number of steps after which print to logfile
335 //! Number of steps after which X is output
337 //! Number of steps after which V is output
339 //! Number of steps after which F is output
341 //! Number of steps after which energies printed
343 //! Number of steps after which compressed trj (.xtc,.tng) is output
344 int nstxout_compressed;
345 //! Initial time (ps)
349 //! Whether we use multiple time stepping
351 //! The multiple time stepping levels
352 std::vector<gmx::MtsLevel> mtsLevels;
353 //! Precision of x in compressed trajectory file
354 real x_compression_precision;
355 //! Requested fourier_spacing, when nk? not set
356 real fourier_spacing;
357 //! Number of k vectors in x dimension for fourier methods for long range electrost.
359 //! Number of k vectors in y dimension for fourier methods for long range electrost.
361 //! Number of k vectors in z dimension for fourier methods for long range electrost.
363 //! Interpolation order for PME
365 //! Real space tolerance for Ewald, determines the real/reciprocal space relative weight
367 //! Real space tolerance for LJ-Ewald
369 //! Normal/3D ewald, or pseudo-2D LR corrections
371 //! Epsilon for PME dipole correction
372 real epsilon_surface;
373 //! Type of combination rule in LJ-PME
374 int ljpme_combination_rule;
375 //! Type of periodic boundary conditions
377 //! Periodic molecules
379 //! Continuation run: starting state is correct (ie. constrained)
380 gmx_bool bContinuation;
381 //! Temperature coupling
383 //! Interval in steps for temperature coupling
385 //! Whether to print nose-hoover chains
386 gmx_bool bPrintNHChains;
387 //! Pressure coupling
389 //! Pressure coupling type
391 //! Interval in steps for pressure coupling
393 //! Pressure coupling time (ps)
395 //! Reference pressure (kJ/(mol nm^3))
397 //! Compressibility ((mol nm^3)/kJ)
399 //! How to scale absolute reference coordinates
400 int refcoord_scaling;
401 //! The COM of the posres atoms
403 //! The B-state COM of the posres atoms
405 //! Random seed for Andersen thermostat (obsolete)
407 //! Per atom pair energy drift tolerance (kJ/mol/ps/atom) for list buffer
409 //! Short range pairlist cut-off (nm)
411 //! Radius for test particle insertion
413 //! Type of electrostatics treatment
415 //! Modify the Coulomb interaction
416 int coulomb_modifier;
417 //! Coulomb switch range start (nm)
418 real rcoulomb_switch;
419 //! Coulomb cutoff (nm)
421 //! Relative dielectric constant
423 //! Relative dielectric constant of the RF
425 //! Always false (no longer supported)
426 bool implicit_solvent;
427 //! Type of Van der Waals treatment
429 //! Modify the Van der Waals interaction
431 //! Van der Waals switch range start (nm)
433 //! Van der Waals cutoff (nm)
435 //! Perform Long range dispersion corrections
437 //! Extension of the table beyond the cut-off, as well as the table length for 1-4 interac.
439 //! Tolerance for shake
441 //! Free energy calculations
443 //! Data for the FEP state
445 //! Whether to do simulated tempering
447 //! Variables for simulated tempering
448 t_simtemp* simtempvals;
449 //! Whether expanded ensembles are used
451 //! Expanded ensemble parameters
452 t_expanded* expandedvals;
453 //! Type of distance restraining
455 //! Force constant for time averaged distance restraints
457 //! Type of weighting of pairs in one restraints
459 //! Use combination of time averaged and instantaneous violations
460 gmx_bool bDisreMixed;
461 //! Frequency of writing pair distances to enx
463 //! Time constant for memory function in disres
465 //! Force constant for orientational restraints
467 //! Time constant for memory function in orires
469 //! Frequency of writing tr(SD) to energy output
471 //! The stepsize for updating
475 //! Number of iterations for convergence of steepest descent in relax_shells
477 //! Stepsize for directional minimization in relax_shells
479 //! Number of steps after which a steepest descents step is done while doing cg
481 //! Number of corrections to the Hessian to keep
483 //! Type of constraint algorithm
485 //! Order of the LINCS Projection Algorithm
487 //! Warn if any bond rotates more than this many degrees
489 //! Number of iterations in the final LINCS step
491 //! Use successive overrelaxation for shake
493 //! Friction coefficient for BD (amu/ps)
495 //! Random seed for SD and BD
497 //! The number of walls
499 //! The type of walls
501 //! The potentail is linear for r<=wall_r_linpot
503 //! The atom type for walls
504 int wall_atomtype[2];
505 //! Number density for walls
506 real wall_density[2];
507 //! Scaling factor for the box for Ewald
508 real wall_ewald_zfac;
510 /* COM pulling data */
511 //! Do we do COM pulling?
513 //! The data for center of mass pulling
514 std::unique_ptr<pull_params_t> pull;
517 //! Whether to use AWH biasing for PMF calculations
519 //! AWH biasing parameters
520 gmx::AwhParams* awhParams;
522 /* Enforced rotation data */
523 //! Whether to calculate enforced rotation potential(s)
525 //! The data for enforced rotation potentials
528 //! Whether to do ion/water position exchanges (CompEL)
530 //! Swap data structure.
533 //! Whether the tpr makes an interactive MD session possible.
535 //! Interactive molecular dynamics
538 //! Acceleration for viscosity calculation
540 //! Triclinic deformation velocities (nm/ps)
542 /*! \brief User determined parameters */
555 //! QM/MM calculation
558 /* Fields for removed features go here (better caching) */
559 //! Whether AdResS is enabled - always false if a valid .tpr was read
561 //! Whether twin-range scheme is active - always false if a valid .tpr was read
562 gmx_bool useTwinRange;
563 //! Whether we have constant acceleration - removed in GROMACS 2022
564 bool useConstantAcceleration;
566 //! KVT object that contains input parameters converted to the new style.
567 gmx::KeyValueTreeObject* params;
569 //! KVT for storing simulation parameters that are not part of the mdp file.
570 std::unique_ptr<gmx::KeyValueTreeObject> internalParameters;
573 int ir_optimal_nstcalcenergy(const t_inputrec* ir);
575 int tcouple_min_integration_steps(int etc);
577 int ir_optimal_nsttcouple(const t_inputrec* ir);
579 int pcouple_min_integration_steps(int epc);
581 int ir_optimal_nstpcouple(const t_inputrec* ir);
583 /* Returns if the Coulomb force or potential is switched to zero */
584 gmx_bool ir_coulomb_switched(const t_inputrec* ir);
586 /* Returns if the Coulomb interactions are zero beyond the rcoulomb.
587 * Note: always returns TRUE for the Verlet cut-off scheme.
589 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir);
591 /* As ir_coulomb_is_zero_at_cutoff, but also returns TRUE for user tabulated
592 * interactions, since these might be zero beyond rcoulomb.
594 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir);
596 /* Returns if the Van der Waals force or potential is switched to zero */
597 gmx_bool ir_vdw_switched(const t_inputrec* ir);
599 /* Returns if the Van der Waals interactions are zero beyond the rvdw.
600 * Note: always returns TRUE for the Verlet cut-off scheme.
602 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir);
604 /* As ir_vdw_is_zero_at_cutoff, but also returns TRUE for user tabulated
605 * interactions, since these might be zero beyond rvdw.
607 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir);
609 /*! \brief Free memory from input record.
611 * All arrays and pointers will be freed.
613 * \param[in] ir The data structure
615 void done_inputrec(t_inputrec* ir);
617 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, gmx_bool bMDPformat);
619 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol);
621 void comp_pull_AB(FILE* fp, const pull_params_t& pull, real ftol, real abstol);
624 gmx_bool inputrecDeform(const t_inputrec* ir);
626 gmx_bool inputrecDynamicBox(const t_inputrec* ir);
628 gmx_bool inputrecPreserveShape(const t_inputrec* ir);
630 gmx_bool inputrecNeedMutot(const t_inputrec* ir);
632 gmx_bool inputrecTwinRange(const t_inputrec* ir);
634 gmx_bool inputrecExclForces(const t_inputrec* ir);
636 gmx_bool inputrecNptTrotter(const t_inputrec* ir);
638 gmx_bool inputrecNvtTrotter(const t_inputrec* ir);
640 gmx_bool inputrecNphTrotter(const t_inputrec* ir);
642 /*! \brief Return true if the simulation is 2D periodic with two walls. */
643 bool inputrecPbcXY2Walls(const t_inputrec* ir);
645 /*! \brief Returns true for MD integator with T and/or P-coupling that supports
646 * calculating a conserved energy quantity.
648 * Note that dynamical integrators without T and P coupling (ie NVE)
649 * return false, i.e. the return value refers to whether there
650 * is a conserved quantity different than the total energy.
652 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir);
654 /*! \brief Returns true when temperature is coupled or constant. */
655 bool integratorHasReferenceTemperature(const t_inputrec* ir);
657 /*! \brief Return the number of bounded dimensions
659 * \param[in] ir The input record with MD parameters
660 * \return the number of dimensions in which
661 * the coordinates of the particles are bounded, starting at X.
663 int inputrec2nboundeddim(const t_inputrec* ir);
665 /*! \brief Returns the number of degrees of freedom in center of mass motion
667 * \param[in] ir The inputrec structure
668 * \return the number of degrees of freedom of the center of mass
670 int ndof_com(const t_inputrec* ir);
672 /*! \brief Returns the maximum reference temperature over all coupled groups
674 * Returns 0 for energy minimization and normal mode computation.
675 * Returns -1 for MD without temperature coupling.
677 * \param[in] ir The inputrec structure
679 real maxReferenceTemperature(const t_inputrec& ir);
681 /*! \brief Returns whether there is an Ewald surface contribution
683 bool haveEwaldSurfaceContribution(const t_inputrec& ir);
685 /*! \brief Check if calculation of the specific FEP type was requested.
687 * \param[in] ir Input record.
688 * \param[in] fepType Free-energy perturbation type to check for.
690 * \returns If the \p fepType is perturbed in this run.
692 bool haveFreeEnergyType(const t_inputrec& ir, int fepType);
694 #endif /* GMX_MDTYPES_INPUTREC_H */