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