Add gapsys softcore function.
[alexxy/gromacs.git] / src / gromacs / mdtypes / inputrec.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #ifndef GMX_MDTYPES_INPUTREC_H
39 #define GMX_MDTYPES_INPUTREC_H
40
41 #include <cstdio>
42
43 #include <memory>
44 #include <vector>
45
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"
50
51 #define EGP_EXCL (1 << 0)
52 #define EGP_TABLE (1 << 1)
53
54 struct gmx_enfrot;
55 struct gmx_enfrotgrp;
56 struct pull_params_t;
57
58 namespace gmx
59 {
60 class Awh;
61 class AwhParams;
62 class KeyValueTreeObject;
63 struct MtsLevel;
64 } // namespace gmx
65
66 struct t_grpopts
67 {
68     //! Number of T-Coupl groups
69     int ngtc = 0;
70     //! Number of of Nose-Hoover chains per group
71     int nhchainlength = 0;
72     //! Number of Freeze groups
73     int ngfrz = 0;
74     //! Number of Energy groups
75     int ngener = 0;
76     //! Number of degrees of freedom in a group
77     real* nrdf = nullptr;
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;
88     //! Tau coupling time
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;
94
95     /* QMMM stuff */
96     //! Number of QM groups
97     int ngQM = 0;
98 };
99
100 struct t_simtemp
101 {
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;
110 };
111
112 struct t_lambda
113 {
114     //! The frequency for calculating dhdl
115     int nstdhdl = 0;
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
125     int n_lambda = 0;
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
135     real sc_alpha = 0;
136     //! Lambda power for soft-core interactions
137     int sc_power = 0;
138     //! R power for soft-core interactions
139     real sc_r_power = 0;
140     //! Free energy soft-core sigma when c6 or c12=0
141     real sc_sigma = 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     //! The specific soft-core function to use
147     SoftcoreType softcoreFunction = SoftcoreType::Beutler;
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::EnumerationArray<FreeEnergyPerturbationCouplingType, bool> separate_dvdl;
150     //! Whether to write a separate dhdl.xvg file note: NOT a gmx_bool, but an enum
151     SeparateDhdlFile separate_dhdl_file = SeparateDhdlFile::Default;
152     //! Whether to calculate+write dhdl derivatives note: NOT a gmx_bool, but an enum
153     DhDlDerivativeCalculation dhdl_derivatives = DhDlDerivativeCalculation::Default;
154     //! The maximum table size for the dH histogram
155     int dh_hist_size = 0;
156     //! The spacing for the dH histogram
157     double dh_hist_spacing = 0;
158 };
159
160 struct t_expanded
161 {
162     //! The frequency of expanded ensemble state changes
163     int nstexpanded = 0;
164     //! Which type of move updating do we use for lambda monte carlo (or no for none)
165     LambdaWeightCalculation elamstats = LambdaWeightCalculation::Default;
166     //! What move set will be we using for state space moves
167     LambdaMoveCalculation elmcmove = LambdaMoveCalculation::Default;
168     //! The method we use to decide of we have equilibrated the weights
169     LambdaWeightWillReachEquilibrium elmceq = LambdaWeightWillReachEquilibrium::Default;
170     //! The minumum number of samples at each lambda for deciding whether we have reached a minimum
171     int equil_n_at_lam = 0;
172     //! Wang-Landau delta at which we stop equilibrating weights
173     real equil_wl_delta = 0;
174     //! Use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating
175     real equil_ratio = 0;
176     //! After equil_steps steps we stop equilibrating the weights
177     int equil_steps = 0;
178     //! After equil_samples total samples (steps/nstfep), we stop equilibrating the weights
179     int equil_samples = 0;
180     //! Random number seed for lambda mc switches
181     int lmc_seed = 0;
182     //! Whether to use minumum variance weighting
183     bool minvar = false;
184     //! The number of samples needed before kicking into minvar routine
185     int minvarmin = 0;
186     //! The offset for the variance in MinVar
187     real minvar_const = 0;
188     //! Range of cvalues used for BAR
189     int c_range = 0;
190     //! Whether to print symmetrized matrices
191     bool bSymmetrizedTMatrix = false;
192     //! How frequently to print the transition matrices
193     int nstTij = 0;
194     //! Number of repetitions in the MC lambda jumps MRS -- VERIFY THIS
195     int lmc_repeats = 0;
196     //! Minimum number of samples for each state before free sampling MRS -- VERIFY THIS!
197     int lmc_forced_nstart = 0;
198     //! Distance in lambda space for the gibbs interval
199     int gibbsdeltalam = 0;
200     //! Scaling factor for Wang-Landau
201     real wl_scale = 0;
202     //! Ratio between largest and smallest number for freezing the weights
203     real wl_ratio = 0;
204     //! Starting delta for Wang-Landau
205     real init_wl_delta = 0;
206     //! Use one over t convergence for Wang-Landau when the delta get sufficiently small
207     bool bWLoneovert = false;
208     //! Did we initialize the weights? TODO: REMOVE FOR 5.0, no longer needed with new logic
209     bool bInit_weights = false;
210     //! To override the main temperature, or define it if it's not defined
211     real mc_temp = 0;
212     //! User-specified initial weights to start with
213     std::vector<real> init_lambda_weights;
214 };
215
216 struct t_rotgrp
217 {
218     //! Rotation type for this group
219     EnforcedRotationGroupType eType;
220     //! Use mass-weighed positions?
221     bool bMassW;
222     //! Number of atoms in the group
223     int nat;
224     //! The global atoms numbers
225     int* ind;
226     //! The reference positions
227     rvec* x_ref;
228     //! The normalized rotation vector
229     rvec inputVec;
230     //! Rate of rotation (degree/ps)
231     real rate;
232     //! Force constant (kJ/(mol nm^2)
233     real k;
234     //! Pivot point of rotation axis (nm)
235     rvec pivot;
236     //! Type of fit to determine actual group angle
237     RotationGroupFitting eFittype;
238     //! Number of angles around the reference angle for which the rotation potential is also evaluated (for fit type 'potential' only)
239     int PotAngle_nstep;
240     //! Distance between two angles in degrees (for fit type 'potential' only)
241     real PotAngle_step;
242     //! Slab distance (nm)
243     real slab_dist;
244     //! Minimum value the gaussian must have so that the force is actually evaluated
245     real min_gaussian;
246     //! Additive constant for radial motion2 and flexible2 potentials (nm^2)
247     real eps;
248 };
249
250 struct t_rot
251 {
252     //! Number of rotation groups
253     int ngrp;
254     //! Output frequency for main rotation outfile
255     int nstrout;
256     //! Output frequency for per-slab data
257     int nstsout;
258     //! Groups to rotate
259     t_rotgrp* grp;
260 };
261
262 struct t_IMD
263 {
264     //! Number of interactive atoms
265     int nat;
266     //! The global indices of the interactive atoms
267     int* ind;
268 };
269
270 struct t_swapGroup
271 {
272     //! Name of the swap group, e.g. NA, CL, SOL
273     char* molname;
274     //! Number of atoms in this group
275     int nat;
276     //! The global ion group atoms numbers
277     int* ind;
278     //! Requested number of molecules of this type per compartment
279     gmx::EnumerationArray<Compartment, int> nmolReq;
280 };
281
282 struct t_swapcoords
283 {
284     //! Period between when a swap is attempted
285     int nstswap;
286     //! Use mass-weighted positions in split group
287     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 */
291     /**@{*/
292     real cyl0r, cyl1r;
293     real cyl0u, cyl1u;
294     real cyl0l, cyl1l;
295     /**@}*/
296     //! Coupling constant (number of swap attempt steps)
297     int nAverage;
298     //! Ion counts may deviate from the requested values by +-threshold before a swap is done
299     real threshold;
300     //! Offset of the swap layer (='bulk') with respect to the compartment-defining layers
301     gmx::EnumerationArray<Compartment, real> bulkOffset;
302     //! Number of groups to be controlled
303     int ngrp;
304     //! All swap groups, including split and solvent
305     t_swapGroup* grp;
306 };
307
308 struct t_inputrec // NOLINT (clang-analyzer-optin.performance.Padding)
309 {
310     t_inputrec();
311     explicit t_inputrec(const t_inputrec&) = delete;
312     t_inputrec& operator=(const t_inputrec&) = delete;
313     ~t_inputrec();
314
315     //! Integration method
316     IntegrationAlgorithm eI = IntegrationAlgorithm::Default;
317     //! Number of steps to be taken
318     int64_t nsteps = 0;
319     //! Used in checkpointing to separate chunks
320     int simulation_part = 0;
321     //! Start at a stepcount >0 (used w. convert-tpr)
322     int64_t init_step = 0;
323     //! Frequency of energy calc. and T/P coupl. upd.
324     int nstcalcenergy = 0;
325     //! Group or verlet cutoffs
326     CutoffScheme cutoff_scheme = CutoffScheme::Default;
327     //! Number of steps before pairlist is generated
328     int nstlist = 0;
329     //! Number of steps after which center of mass motion is removed
330     int nstcomm = 0;
331     //! Center of mass motion removal algorithm
332     ComRemovalAlgorithm comm_mode = ComRemovalAlgorithm::Default;
333     //! Number of steps after which print to logfile
334     int nstlog = 0;
335     //! Number of steps after which X is output
336     int nstxout = 0;
337     //! Number of steps after which V is output
338     int nstvout = 0;
339     //! Number of steps after which F is output
340     int nstfout = 0;
341     //! Number of steps after which energies printed
342     int nstenergy = 0;
343     //! Number of steps after which compressed trj (.xtc,.tng) is output
344     int nstxout_compressed = 0;
345     //! Initial time (ps)
346     double init_t = 0;
347     //! Time step (ps)
348     double delta_t = 0;
349     //! Whether we use multiple time stepping
350     bool useMts = false;
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 = 0;
355     //! Requested fourier_spacing, when nk? not set
356     real fourier_spacing = 0;
357     //! Number of k vectors in x dimension for fourier methods for long range electrost.
358     int nkx = 0;
359     //! Number of k vectors in y dimension for fourier methods for long range electrost.
360     int nky = 0;
361     //! Number of k vectors in z dimension for fourier methods for long range electrost.
362     int nkz = 0;
363     //! Interpolation order for PME
364     int pme_order = 0;
365     //! Real space tolerance for Ewald, determines the real/reciprocal space relative weight
366     real ewald_rtol = 0;
367     //! Real space tolerance for LJ-Ewald
368     real ewald_rtol_lj = 0;
369     //! Normal/3D ewald, or pseudo-2D LR corrections
370     EwaldGeometry ewald_geometry = EwaldGeometry::Default;
371     //! Epsilon for PME dipole correction
372     real epsilon_surface = 0;
373     //! Type of combination rule in LJ-PME
374     LongRangeVdW ljpme_combination_rule = LongRangeVdW::Default;
375     //! Type of periodic boundary conditions
376     PbcType pbcType = PbcType::Default;
377     //! Periodic molecules
378     bool bPeriodicMols = false;
379     //! Continuation run: starting state is correct (ie. constrained)
380     bool bContinuation = false;
381     //! Temperature coupling
382     TemperatureCoupling etc = TemperatureCoupling::Default;
383     //! Interval in steps for temperature coupling
384     int nsttcouple = 0;
385     //! Whether to print nose-hoover chains
386     bool bPrintNHChains = false;
387     //! Pressure coupling
388     PressureCoupling epc = PressureCoupling::Default;
389     //! Pressure coupling type
390     PressureCouplingType epct = PressureCouplingType::Default;
391     //! Interval in steps for pressure coupling
392     int nstpcouple = 0;
393     //! Pressure coupling time (ps)
394     real tau_p = 0;
395     //! Reference pressure (kJ/(mol nm^3))
396     tensor ref_p = { { 0 } };
397     //! Compressibility ((mol nm^3)/kJ)
398     tensor compress = { { 0 } };
399     //! How to scale absolute reference coordinates
400     RefCoordScaling refcoord_scaling = RefCoordScaling::Default;
401     //! The COM of the posres atoms
402     rvec posres_com = { 0, 0, 0 };
403     //! The B-state COM of the posres atoms
404     rvec posres_comB = { 0, 0, 0 };
405     //! Random seed for Andersen thermostat (obsolete)
406     int andersen_seed = 0;
407     //! Per atom pair energy drift tolerance (kJ/mol/ps/atom) for list buffer
408     real verletbuf_tol = 0;
409     //! Short range pairlist cut-off (nm)
410     real rlist = 0;
411     //! Radius for test particle insertion
412     real rtpi = 0;
413     //! Type of electrostatics treatment
414     CoulombInteractionType coulombtype = CoulombInteractionType::Default;
415     //! Modify the Coulomb interaction
416     InteractionModifiers coulomb_modifier = InteractionModifiers::Default;
417     //! Coulomb switch range start (nm)
418     real rcoulomb_switch = 0;
419     //! Coulomb cutoff (nm)
420     real rcoulomb = 0;
421     //! Relative dielectric constant
422     real epsilon_r = 0;
423     //! Relative dielectric constant of the RF
424     real epsilon_rf = 0;
425     //! Always false (no longer supported)
426     bool implicit_solvent = false;
427     //! Type of Van der Waals treatment
428     VanDerWaalsType vdwtype = VanDerWaalsType::Default;
429     //! Modify the Van der Waals interaction
430     InteractionModifiers vdw_modifier = InteractionModifiers::Default;
431     //! Van der Waals switch range start (nm)
432     real rvdw_switch = 0;
433     //! Van der Waals cutoff (nm)
434     real rvdw = 0;
435     //! Perform Long range dispersion corrections
436     DispersionCorrectionType eDispCorr = DispersionCorrectionType::Default;
437     //! Extension of the table beyond the cut-off, as well as the table length for 1-4 interac.
438     real tabext = 0;
439     //! Tolerance for shake
440     real shake_tol = 0;
441     //! Free energy calculations
442     FreeEnergyPerturbationType efep = FreeEnergyPerturbationType::Default;
443     //! Data for the FEP state
444     std::unique_ptr<t_lambda> fepvals;
445     //! Whether to do simulated tempering
446     bool bSimTemp = false;
447     //! Variables for simulated tempering
448     std::unique_ptr<t_simtemp> simtempvals;
449     //! Whether expanded ensembles are used
450     bool bExpanded = false;
451     //! Expanded ensemble parameters
452     std::unique_ptr<t_expanded> expandedvals;
453     //! Type of distance restraining
454     DistanceRestraintRefinement eDisre = DistanceRestraintRefinement::Default;
455     //! Force constant for time averaged distance restraints
456     real dr_fc = 0;
457     //! Type of weighting of pairs in one restraints
458     DistanceRestraintWeighting eDisreWeighting = DistanceRestraintWeighting::Default;
459     //! Use combination of time averaged and instantaneous violations
460     bool bDisreMixed = false;
461     //! Frequency of writing pair distances to enx
462     int nstdisreout = 0;
463     //! Time constant for memory function in disres
464     real dr_tau = 0;
465     //! Force constant for orientational restraints
466     real orires_fc = 0;
467     //! Time constant for memory function in orires
468     real orires_tau = 0;
469     //! Frequency of writing tr(SD) to energy output
470     int nstorireout = 0;
471     //! The stepsize for updating
472     real em_stepsize = 0;
473     //! The tolerance
474     real em_tol = 0;
475     //! Number of iterations for convergence of steepest descent in relax_shells
476     int niter = 0;
477     //! Stepsize for directional minimization in relax_shells
478     real fc_stepsize = 0;
479     //! Number of steps after which a steepest descents step is done while doing cg
480     int nstcgsteep = 0;
481     //! Number of corrections to the Hessian to keep
482     int nbfgscorr = 0;
483     //! Type of constraint algorithm
484     ConstraintAlgorithm eConstrAlg = ConstraintAlgorithm::Default;
485     //! Order of the LINCS Projection Algorithm
486     int nProjOrder = 0;
487     //! Warn if any bond rotates more than this many degrees
488     real LincsWarnAngle = 0;
489     //! Number of iterations in the final LINCS step
490     int nLincsIter = 0;
491     //! Use successive overrelaxation for shake
492     bool bShakeSOR = false;
493     //! Friction coefficient for BD (amu/ps)
494     real bd_fric = 0;
495     //! Random seed for SD and BD
496     int64_t ld_seed = 0;
497     //! The number of walls
498     int nwall = 0;
499     //! The type of walls
500     WallType wall_type = WallType::Default;
501     //! The potentail is linear for r<=wall_r_linpot
502     real wall_r_linpot = 0;
503     //! The atom type for walls
504     int wall_atomtype[2] = { 0, 0 };
505     //! Number density for walls
506     real wall_density[2] = { 0, 0 };
507     //! Scaling factor for the box for Ewald
508     real wall_ewald_zfac = 0;
509
510     /* COM pulling data */
511     //! Do we do COM pulling?
512     bool bPull = false;
513     //! The data for center of mass pulling
514     std::unique_ptr<pull_params_t> pull;
515
516     /* AWH bias data */
517     //! Whether to use AWH biasing for PMF calculations
518     bool bDoAwh = false;
519     //! AWH biasing parameters
520     std::unique_ptr<gmx::AwhParams> awhParams;
521
522     /* Enforced rotation data */
523     //! Whether to calculate enforced rotation potential(s)
524     bool bRot = false;
525     //! The data for enforced rotation potentials
526     t_rot* rot = nullptr;
527
528     //! Whether to do ion/water position exchanges (CompEL)
529     SwapType eSwapCoords = SwapType::Default;
530     //! Swap data structure.
531     t_swapcoords* swap = nullptr;
532
533     //! Whether the tpr makes an interactive MD session possible.
534     bool bIMD = false;
535     //! Interactive molecular dynamics
536     t_IMD* imd = nullptr;
537
538     //! Acceleration for viscosity calculation
539     real cos_accel = 0;
540     //! Triclinic deformation velocities (nm/ps)
541     tensor deform = { { 0 } };
542     /*! \brief User determined parameters */
543     /**@{*/
544     int  userint1  = 0;
545     int  userint2  = 0;
546     int  userint3  = 0;
547     int  userint4  = 0;
548     real userreal1 = 0;
549     real userreal2 = 0;
550     real userreal3 = 0;
551     real userreal4 = 0;
552     /**@}*/
553     //! Group options
554     t_grpopts opts;
555     //! QM/MM calculation
556     bool bQMMM = false;
557
558     /* Fields for removed features go here (better caching) */
559     //! Whether AdResS is enabled - always false if a valid .tpr was read
560     bool bAdress = false;
561     //! Whether twin-range scheme is active - always false if a valid .tpr was read
562     bool useTwinRange = false;
563     //! Whether we have constant acceleration - removed in GROMACS 2022
564     bool useConstantAcceleration = false;
565
566     //! KVT object that contains input parameters converted to the new style.
567     gmx::KeyValueTreeObject* params = nullptr;
568
569     //! KVT for storing simulation parameters that are not part of the mdp file.
570     std::unique_ptr<gmx::KeyValueTreeObject> internalParameters;
571 };
572
573 int ir_optimal_nstcalcenergy(const t_inputrec* ir);
574
575 int tcouple_min_integration_steps(TemperatureCoupling etc);
576
577 int ir_optimal_nsttcouple(const t_inputrec* ir);
578
579 int pcouple_min_integration_steps(PressureCoupling epc);
580
581 int ir_optimal_nstpcouple(const t_inputrec* ir);
582
583 /* Returns if the Coulomb force or potential is switched to zero */
584 bool ir_coulomb_switched(const t_inputrec* ir);
585
586 /* Returns if the Coulomb interactions are zero beyond the rcoulomb.
587  * Note: always returns TRUE for the Verlet cut-off scheme.
588  */
589 bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir);
590
591 /* As ir_coulomb_is_zero_at_cutoff, but also returns TRUE for user tabulated
592  * interactions, since these might be zero beyond rcoulomb.
593  */
594 bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir);
595
596 /* Returns if the Van der Waals force or potential is switched to zero */
597 bool ir_vdw_switched(const t_inputrec* ir);
598
599 /* Returns if the Van der Waals interactions are zero beyond the rvdw.
600  * Note: always returns TRUE for the Verlet cut-off scheme.
601  */
602 bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir);
603
604 /* As ir_vdw_is_zero_at_cutoff, but also returns TRUE for user tabulated
605  * interactions, since these might be zero beyond rvdw.
606  */
607 bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir);
608
609 /*! \brief Free memory from input record.
610  *
611  * All arrays and pointers will be freed.
612  *
613  * \param[in] ir The data structure
614  */
615 void done_inputrec(t_inputrec* ir);
616
617 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, bool bMDPformat);
618
619 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol);
620
621 void comp_pull_AB(FILE* fp, const pull_params_t& pull, real ftol, real abstol);
622
623
624 bool inputrecDeform(const t_inputrec* ir);
625
626 bool inputrecDynamicBox(const t_inputrec* ir);
627
628 bool inputrecPreserveShape(const t_inputrec* ir);
629
630 bool inputrecNeedMutot(const t_inputrec* ir);
631
632 bool inputrecTwinRange(const t_inputrec* ir);
633
634 bool inputrecExclForces(const t_inputrec* ir);
635
636 bool inputrecNptTrotter(const t_inputrec* ir);
637
638 bool inputrecNvtTrotter(const t_inputrec* ir);
639
640 bool inputrecNphTrotter(const t_inputrec* ir);
641
642 /*! \brief Return true if the simulation is 2D periodic with two walls. */
643 bool inputrecPbcXY2Walls(const t_inputrec* ir);
644
645 //! \brief Return true if the simulation has frozen atoms (non-trivial freeze groups).
646 bool inputrecFrozenAtoms(const t_inputrec* ir);
647
648 /*! \brief Returns true for MD integator with T and/or P-coupling that supports
649  * calculating a conserved energy quantity.
650  *
651  * Note that dynamical integrators without T and P coupling (ie NVE)
652  * return false, i.e. the return value refers to whether there
653  * is a conserved quantity different than the total energy.
654  */
655 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir);
656
657 /*! \brief Returns true when temperature is coupled or constant. */
658 bool integratorHasReferenceTemperature(const t_inputrec* ir);
659
660 /*! \brief Return the number of bounded dimensions
661  *
662  * \param[in] ir The input record with MD parameters
663  * \return the number of dimensions in which
664  * the coordinates of the particles are bounded, starting at X.
665  */
666 int inputrec2nboundeddim(const t_inputrec* ir);
667
668 /*! \brief Returns the number of degrees of freedom in center of mass motion
669  *
670  * \param[in] ir  The inputrec structure
671  * \return the number of degrees of freedom of the center of mass
672  */
673 int ndof_com(const t_inputrec* ir);
674
675 /*! \brief Returns the maximum reference temperature over all coupled groups
676  *
677  * Returns 0 for energy minimization and normal mode computation.
678  * Returns -1 for MD without temperature coupling.
679  *
680  * \param[in] ir  The inputrec structure
681  */
682 real maxReferenceTemperature(const t_inputrec& ir);
683
684 /*! \brief Returns whether there is an Ewald surface contribution
685  */
686 bool haveEwaldSurfaceContribution(const t_inputrec& ir);
687
688 /*! \brief Check if calculation of the specific FEP type was requested.
689  *
690  * \param[in] ir       Input record.
691  * \param[in] fepType  Free-energy perturbation type to check for.
692  *
693  * \returns If the \p fepType is perturbed in this run.
694  */
695 bool haveFreeEnergyType(const t_inputrec& ir, int fepType);
696
697 #endif /* GMX_MDTYPES_INPUTREC_H */