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,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_MDTYPES_INPUTREC_H
38 #define GMX_MDTYPES_INPUTREC_H
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/mdtypes/md_enums.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/real.h"
49 #define EGP_EXCL (1<<0)
50 #define EGP_TABLE (1<<1)
57 typedef struct t_swap
*gmx_swapcoords_t
;
63 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
102 //! Level of theory in the QM calculation
104 //! Basisset in the QM calculation
106 //! Total charge in the QM region
108 //! Spin multiplicicty in the QM region
110 //! Surface hopping (diabatic hop only)
112 //! Number of orbiatls in the active space
114 //! Number of electrons in the active space
116 //! At which gap (A.U.) the SA is switched on
118 //! At which gap (A.U.) the SA is switched off
120 //! In how many steps SA goes from 1-1 to 0.5-0.5
126 //! Simulated temperature scaling; linear or exponential
128 //! The low temperature for simulated tempering
130 //! The high temperature for simulated tempering
132 //! The range of temperatures used for simulated tempering
138 //! The frequency for calculating dhdl
140 //! 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)
142 //! The initial number of the state
144 //! Change of lambda per time step (fraction of (0.1)
146 //! Print no, total or potential energies in dhdl
147 int edHdLPrintEnergy
;
148 //! The number of foreign lambda points
150 //! The array of all lambda values
152 //! The number of neighboring lambda states to calculate the energy for in up and down directions (-1 for all)
153 int lambda_neighbors
;
154 //! The first lambda to calculate energies for
156 //! The last lambda +1 to calculate energies for
158 //! Free energy soft-core parameter
160 //! Lambda power for soft-core interactions
162 //! R power for soft-core interactions
164 //! Free energy soft-core sigma when c6 or c12=0
166 //! Free energy soft-core sigma for ?????
168 //! Use softcore for the coulomb portion as well (default FALSE)
170 //! Whether to print the dvdl term associated with this term; if it is not specified as separate, it is lumped with the FEP term
171 gmx_bool separate_dvdl
[efptNR
];
172 //! Whether to write a separate dhdl.xvg file note: NOT a gmx_bool, but an enum
173 int separate_dhdl_file
;
174 //! Whether to calculate+write dhdl derivatives note: NOT a gmx_bool, but an enum
175 int dhdl_derivatives
;
176 //! The maximum table size for the dH histogram
178 //! The spacing for the dH histogram
179 double dh_hist_spacing
;
183 //! The frequency of expanded ensemble state changes
185 //! Which type of move updating do we use for lambda monte carlo (or no for none)
187 //! What move set will be we using for state space moves
189 //! The method we use to decide of we have equilibrated the weights
191 //! The minumum number of samples at each lambda for deciding whether we have reached a minimum
193 //! Wang-Landau delta at which we stop equilibrating weights
195 //! Use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating
197 //! After equil_steps steps we stop equilibrating the weights
199 //! After equil_samples total samples (steps/nstfep), we stop equilibrating the weights
201 //! Random number seed for lambda mc switches
203 //! Whether to use minumum variance weighting
205 //! The number of samples needed before kicking into minvar routine
207 //! The offset for the variance in MinVar
209 //! Range of cvalues used for BAR
211 //! Whether to print symmetrized matrices
212 gmx_bool bSymmetrizedTMatrix
;
213 //! How frequently to print the transition matrices
215 //! Number of repetitions in the MC lambda jumps MRS -- VERIFY THIS
217 //! Minimum number of samples for each state before free sampling MRS -- VERIFY THIS!
218 int lmc_forced_nstart
;
219 //! Distance in lambda space for the gibbs interval
221 //! Scaling factor for Wang-Landau
223 //! Ratio between largest and smallest number for freezing the weights
225 //! Starting delta for Wang-Landau
227 //! Use one over t convergence for Wang-Landau when the delta get sufficiently small
228 gmx_bool bWLoneovert
;
229 //! Did we initialize the weights? TODO: REMOVE FOR 5.0, no longer needed with new logic
230 gmx_bool bInit_weights
;
231 //! To override the main temperature, or define it if it's not defined
233 //! User-specified initial weights to start with
234 real
*init_lambda_weights
;
239 //! Rotation type for this group
241 //! Use mass-weighed positions?
243 //! Number of atoms in the group
245 //! The global atoms numbers
247 //! The reference positions
249 //! The normalized rotation vector
251 //! Rate of rotation (degree/ps)
253 //! Force constant (kJ/(mol nm^2)
255 //! Pivot point of rotation axis (nm)
257 //! Type of fit to determine actual group angle
259 //! Number of angles around the reference angle for which the rotation potential is also evaluated (for fit type 'potential' only)
261 //! Distance between two angles in degrees (for fit type 'potential' only)
263 //! Slab distance (nm)
265 //! Minimum value the gaussian must have so that the force is actually evaluated
267 //! Additive constant for radial motion2 and flexible2 potentials (nm^2)
273 //! Number of rotation groups
275 //! Output frequency for main rotation outfile
277 //! Output frequency for per-slab data
285 //! Number of interactive atoms
287 //! The global indices of the interactive atoms
289 //! Stores non-inputrec IMD data
295 //! Name of the swap group, e.g. NA, CL, SOL
297 //! Number of atoms in this group
299 //! The global ion group atoms numbers
301 //! Requested number of molecules of this type per compartment
302 int nmolReq
[eCompNR
];
307 //! Period between when a swap is attempted
309 //! Use mass-weighted positions in split group
310 gmx_bool massw_split
[2];
311 /*! \brief Split cylinders defined by radius, upper and lower
312 * extension. The split cylinders define the channels and are
313 * each anchored in the center of the split group */
319 //! Coupling constant (number of swap attempt steps)
321 //! Ion counts may deviate from the requested values by +-threshold before a swap is done
323 //! Offset of the swap layer (='bulk') with respect to the compartment-defining layers
324 real bulkOffset
[eCompNR
];
325 //! Number of groups to be controlled
327 //! All swap groups, including split and solvent
329 //! Swap private data accessible in swapcoords.cpp
330 gmx_swapcoords_t si_priv
;
336 explicit t_inputrec(const t_inputrec
&) = delete;
337 t_inputrec
&operator=(const t_inputrec
&) = delete;
340 //! Integration method
342 //! Number of steps to be taken
344 //! Used in checkpointing to separate chunks
346 //! Start at a stepcount >0 (used w. convert-tpr)
348 //! Frequency of energy calc. and T/P coupl. upd.
350 //! Group or verlet cutoffs
352 //! Which neighbor search method should we use?
354 //! Number of steps before pairlist is generated
356 //! Number of cells per rlong
358 //! Number of steps after which center of mass motion is removed
360 //! Center of mass motion removal algorithm
362 //! Number of steps after which print to logfile
364 //! Number of steps after which X is output
366 //! Number of steps after which V is output
368 //! Number of steps after which F is output
370 //! Number of steps after which energies printed
372 //! Number of steps after which compressed trj (.xtc,.tng) is output
373 int nstxout_compressed
;
374 //! Initial time (ps)
378 //! Precision of x in compressed trajectory file
379 real x_compression_precision
;
380 //! Requested fourier_spacing, when nk? not set
381 real fourier_spacing
;
382 //! Number of k vectors in x dimension for fourier methods for long range electrost.
384 //! Number of k vectors in y dimension for fourier methods for long range electrost.
386 //! Number of k vectors in z dimension for fourier methods for long range electrost.
388 //! Interpolation order for PME
390 //! Real space tolerance for Ewald, determines the real/reciprocal space relative weight
392 //! Real space tolerance for LJ-Ewald
394 //! Normal/3D ewald, or pseudo-2D LR corrections
396 //! Epsilon for PME dipole correction
397 real epsilon_surface
;
398 //! Type of combination rule in LJ-PME
399 int ljpme_combination_rule
;
400 //! Type of periodic boundary conditions
402 //! Periodic molecules
404 //! Continuation run: starting state is correct (ie. constrained)
405 gmx_bool bContinuation
;
406 //! Temperature coupling
408 //! Interval in steps for temperature coupling
410 //! Whether to print nose-hoover chains
411 gmx_bool bPrintNHChains
;
412 //! Pressure coupling
414 //! Pressure coupling type
416 //! Interval in steps for pressure coupling
418 //! Pressure coupling time (ps)
420 //! Reference pressure (kJ/(mol nm^3))
422 //! Compressibility ((mol nm^3)/kJ)
424 //! How to scale absolute reference coordinates
425 int refcoord_scaling
;
426 //! The COM of the posres atoms
428 //! The B-state COM of the posres atoms
430 //! Random seed for Andersen thermostat (obsolete)
432 //! Per atom pair energy drift tolerance (kJ/mol/ps/atom) for list buffer
434 //! Short range pairlist cut-off (nm)
436 //! Radius for test particle insertion
438 //! Type of electrostatics treatment
440 //! Modify the Coulomb interaction
441 int coulomb_modifier
;
442 //! Coulomb switch range start (nm)
443 real rcoulomb_switch
;
444 //! Coulomb cutoff (nm)
446 //! Relative dielectric constant
448 //! Relative dielectric constant of the RF
450 //! Always false (no longer supported)
451 bool implicit_solvent
;
452 //! Type of Van der Waals treatment
454 //! Modify the Van der Waals interaction
456 //! Van der Waals switch range start (nm)
458 //! Van der Waals cutoff (nm)
460 //! Perform Long range dispersion corrections
462 //! Extension of the table beyond the cut-off, as well as the table length for 1-4 interac.
464 //! Tolerance for shake
466 //! Free energy calculations
468 //! Data for the FEP state
470 //! Whether to do simulated tempering
472 //! Variables for simulated tempering
473 t_simtemp
*simtempvals
;
474 //! Whether expanded ensembles are used
476 //! Expanded ensemble parameters
477 t_expanded
*expandedvals
;
478 //! Type of distance restraining
480 //! Force constant for time averaged distance restraints
482 //! Type of weighting of pairs in one restraints
484 //! Use combination of time averaged and instantaneous violations
485 gmx_bool bDisreMixed
;
486 //! Frequency of writing pair distances to enx
488 //! Time constant for memory function in disres
490 //! Force constant for orientational restraints
492 //! Time constant for memory function in orires
494 //! Frequency of writing tr(SD) to energy output
496 //! The stepsize for updating
500 //! Number of iterations for convergence of steepest descent in relax_shells
502 //! Stepsize for directional minimization in relax_shells
504 //! Number of steps after which a steepest descents step is done while doing cg
506 //! Number of corrections to the Hessian to keep
508 //! Type of constraint algorithm
510 //! Order of the LINCS Projection Algorithm
512 //! Warn if any bond rotates more than this many degrees
514 //! Number of iterations in the final LINCS step
516 //! Use successive overrelaxation for shake
518 //! Friction coefficient for BD (amu/ps)
520 //! Random seed for SD and BD
522 //! The number of walls
524 //! The type of walls
526 //! The potentail is linear for r<=wall_r_linpot
528 //! The atom type for walls
529 int wall_atomtype
[2];
530 //! Number density for walls
531 real wall_density
[2];
532 //! Scaling factor for the box for Ewald
533 real wall_ewald_zfac
;
535 /* COM pulling data */
536 //! Do we do COM pulling?
538 //! The data for center of mass pulling
540 // TODO: Remove this by converting pull into a ForceProvider
541 //! The COM pull force calculation data structure
545 //! Whether to use AWH biasing for PMF calculations
547 //! AWH biasing parameters
548 gmx::AwhParams
*awhParams
;
550 /* Enforced rotation data */
551 //! Whether to calculate enforced rotation potential(s)
553 //! The data for enforced rotation potentials
556 //! Whether to do ion/water position exchanges (CompEL)
558 //! Swap data structure.
561 //! Whether to do an interactive MD session
563 //! Interactive molecular dynamics
566 //! Acceleration for viscosity calculation
568 //! Triclinic deformation velocities (nm/ps)
570 /*! \brief User determined parameters */
583 //! QM/MM calculation
585 //! Constraints on QM bonds
587 //! Scheme: ONIOM or normal
589 //! Factor for scaling the MM charges in QM calc.
592 /* Fields for removed features go here (better caching) */
593 //! Whether AdResS is enabled - always false if a valid .tpr was read
595 //! Whether twin-range scheme is active - always false if a valid .tpr was read
596 gmx_bool useTwinRange
;
598 //! KVT object that contains input parameters converted to the new style.
599 gmx::KeyValueTreeObject
*params
;
602 int ir_optimal_nstcalcenergy(const t_inputrec
*ir
);
604 int tcouple_min_integration_steps(int etc
);
606 int ir_optimal_nsttcouple(const t_inputrec
*ir
);
608 int pcouple_min_integration_steps(int epc
);
610 int ir_optimal_nstpcouple(const t_inputrec
*ir
);
612 /* Returns if the Coulomb force or potential is switched to zero */
613 gmx_bool
ir_coulomb_switched(const t_inputrec
*ir
);
615 /* Returns if the Coulomb interactions are zero beyond the rcoulomb.
616 * Note: always returns TRUE for the Verlet cut-off scheme.
618 gmx_bool
ir_coulomb_is_zero_at_cutoff(const t_inputrec
*ir
);
620 /* As ir_coulomb_is_zero_at_cutoff, but also returns TRUE for user tabulated
621 * interactions, since these might be zero beyond rcoulomb.
623 gmx_bool
ir_coulomb_might_be_zero_at_cutoff(const t_inputrec
*ir
);
625 /* Returns if the Van der Waals force or potential is switched to zero */
626 gmx_bool
ir_vdw_switched(const t_inputrec
*ir
);
628 /* Returns if the Van der Waals interactions are zero beyond the rvdw.
629 * Note: always returns TRUE for the Verlet cut-off scheme.
631 gmx_bool
ir_vdw_is_zero_at_cutoff(const t_inputrec
*ir
);
633 /* As ir_vdw_is_zero_at_cutoff, but also returns TRUE for user tabulated
634 * interactions, since these might be zero beyond rvdw.
636 gmx_bool
ir_vdw_might_be_zero_at_cutoff(const t_inputrec
*ir
);
638 /*! \brief Free memory from input record.
640 * All arrays and pointers will be freed.
642 * \param[in] ir The data structure
644 void done_inputrec(t_inputrec
*ir
);
646 void pr_inputrec(FILE *fp
, int indent
, const char *title
, const t_inputrec
*ir
,
647 gmx_bool bMDPformat
);
649 void cmp_inputrec(FILE *fp
, const t_inputrec
*ir1
, const t_inputrec
*ir2
, real ftol
, real abstol
);
651 void comp_pull_AB(FILE *fp
, pull_params_t
*pull
, real ftol
, real abstol
);
654 gmx_bool
inputrecDeform(const t_inputrec
*ir
);
656 gmx_bool
inputrecDynamicBox(const t_inputrec
*ir
);
658 gmx_bool
inputrecPreserveShape(const t_inputrec
*ir
);
660 gmx_bool
inputrecNeedMutot(const t_inputrec
*ir
);
662 gmx_bool
inputrecTwinRange(const t_inputrec
*ir
);
664 gmx_bool
inputrecExclForces(const t_inputrec
*ir
);
666 gmx_bool
inputrecNptTrotter(const t_inputrec
*ir
);
668 gmx_bool
inputrecNvtTrotter(const t_inputrec
*ir
);
670 gmx_bool
inputrecNphTrotter(const t_inputrec
*ir
);
672 /*! \brief Return true if the simulation is 2D periodic with two walls. */
673 bool inputrecPbcXY2Walls(const t_inputrec
*ir
);
675 /*! \brief Returns true for MD integator with T and/or P-coupling that supports
676 * calculating the conserved energy quantity.
678 bool integratorHasConservedEnergyQuantity(const t_inputrec
*ir
);
680 /*! \brief Returns true when temperature is coupled or constant. */
681 bool integratorHasReferenceTemperature(const t_inputrec
*ir
);
683 /*! \brief Return the number of bounded dimensions
685 * \param[in] ir The input record with MD parameters
686 * \return the number of dimensions in which
687 * the coordinates of the particles are bounded, starting at X.
689 int inputrec2nboundeddim(const t_inputrec
*ir
);
691 /*! \brief Returns the number of degrees of freedom in center of mass motion
693 * \param[in] ir the inputrec structure
694 * \return the number of degrees of freedom of the center of mass
696 int ndof_com(const t_inputrec
*ir
);
698 #endif /* GMX_MDTYPES_INPUTREC_H */