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