Merge branch release-2019 into release-2020
[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 } // namespace gmx
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 {
181     //! The frequency of expanded ensemble state changes
182     int nstexpanded;
183     //! Which type of move updating do we use for lambda monte carlo (or no for none)
184     int elamstats;
185     //! What move set will be we using for state space moves
186     int elmcmove;
187     //! The method we use to decide of we have equilibrated the weights
188     int elmceq;
189     //! The minumum number of samples at each lambda for deciding whether we have reached a minimum
190     int equil_n_at_lam;
191     //! Wang-Landau delta at which we stop equilibrating weights
192     real equil_wl_delta;
193     //! Use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating
194     real equil_ratio;
195     //! After equil_steps steps we stop equilibrating the weights
196     int equil_steps;
197     //! After equil_samples total samples (steps/nstfep), we stop equilibrating the weights
198     int equil_samples;
199     //! Random number seed for lambda mc switches
200     int lmc_seed;
201     //! Whether to use minumum variance weighting
202     gmx_bool minvar;
203     //! The number of samples needed before kicking into minvar routine
204     int minvarmin;
205     //! The offset for the variance in MinVar
206     real minvar_const;
207     //! Range of cvalues used for BAR
208     int c_range;
209     //! Whether to print symmetrized matrices
210     gmx_bool bSymmetrizedTMatrix;
211     //! How frequently to print the transition matrices
212     int nstTij;
213     //! Number of repetitions in the MC lambda jumps MRS -- VERIFY THIS
214     int lmc_repeats;
215     //! Minimum number of samples for each state before free sampling MRS -- VERIFY THIS!
216     int lmc_forced_nstart;
217     //! Distance in lambda space for the gibbs interval
218     int gibbsdeltalam;
219     //! Scaling factor for Wang-Landau
220     real wl_scale;
221     //! Ratio between largest and smallest number for freezing the weights
222     real wl_ratio;
223     //! Starting delta for Wang-Landau
224     real init_wl_delta;
225     //! Use one over t convergence for Wang-Landau when the delta get sufficiently small
226     gmx_bool bWLoneovert;
227     //! Did we initialize the weights? TODO: REMOVE FOR 5.0, no longer needed with new logic
228     gmx_bool bInit_weights;
229     //! To override the main temperature, or define it if it's not defined
230     real mc_temp;
231     //! User-specified initial weights to start with
232     real* init_lambda_weights;
233 };
234
235 struct t_rotgrp
236 {
237     //! Rotation type for this group
238     int eType;
239     //! Use mass-weighed positions?
240     int bMassW;
241     //! Number of atoms in the group
242     int nat;
243     //! The global atoms numbers
244     int* ind;
245     //! The reference positions
246     rvec* x_ref;
247     //! The normalized rotation vector
248     rvec inputVec;
249     //! Rate of rotation (degree/ps)
250     real rate;
251     //! Force constant (kJ/(mol nm^2)
252     real k;
253     //! Pivot point of rotation axis (nm)
254     rvec pivot;
255     //! Type of fit to determine actual group angle
256     int eFittype;
257     //! Number of angles around the reference angle for which the rotation potential is also evaluated (for fit type 'potential' only)
258     int PotAngle_nstep;
259     //! Distance between two angles in degrees (for fit type 'potential' only)
260     real PotAngle_step;
261     //! Slab distance (nm)
262     real slab_dist;
263     //! Minimum value the gaussian must have so that the force is actually evaluated
264     real min_gaussian;
265     //! Additive constant for radial motion2 and flexible2 potentials (nm^2)
266     real eps;
267 };
268
269 struct t_rot
270 {
271     //! Number of rotation groups
272     int ngrp;
273     //! Output frequency for main rotation outfile
274     int nstrout;
275     //! Output frequency for per-slab data
276     int nstsout;
277     //! Groups to rotate
278     t_rotgrp* grp;
279 };
280
281 struct t_IMD
282 {
283     //! Number of interactive atoms
284     int nat;
285     //! The global indices of the interactive atoms
286     int* ind;
287 };
288
289 struct t_swapGroup
290 {
291     //! Name of the swap group, e.g. NA, CL, SOL
292     char* molname;
293     //! Number of atoms in this group
294     int nat;
295     //! The global ion group atoms numbers
296     int* ind;
297     //! Requested number of molecules of this type per compartment
298     int nmolReq[eCompNR];
299 };
300
301 struct t_swapcoords
302 {
303     //! Period between when a swap is attempted
304     int nstswap;
305     //! Use mass-weighted positions in split group
306     gmx_bool massw_split[2];
307     /*! \brief Split cylinders defined by radius, upper and lower
308      * extension. The split cylinders define the channels and are
309      * each anchored in the center of the split group */
310     /**@{*/
311     real cyl0r, cyl1r;
312     real cyl0u, cyl1u;
313     real cyl0l, cyl1l;
314     /**@}*/
315     //! Coupling constant (number of swap attempt steps)
316     int nAverage;
317     //! Ion counts may deviate from the requested values by +-threshold before a swap is done
318     real threshold;
319     //! Offset of the swap layer (='bulk') with respect to the compartment-defining layers
320     real bulkOffset[eCompNR];
321     //! Number of groups to be controlled
322     int ngrp;
323     //! All swap groups, including split and solvent
324     t_swapGroup* grp;
325 };
326
327 struct t_inputrec // NOLINT (clang-analyzer-optin.performance.Padding)
328 {
329     t_inputrec();
330     explicit t_inputrec(const t_inputrec&) = delete;
331     t_inputrec& operator=(const t_inputrec&) = delete;
332     ~t_inputrec();
333
334     //! Integration method
335     int eI;
336     //! Number of steps to be taken
337     int64_t nsteps;
338     //! Used in checkpointing to separate chunks
339     int simulation_part;
340     //! Start at a stepcount >0 (used w. convert-tpr)
341     int64_t init_step;
342     //! Frequency of energy calc. and T/P coupl. upd.
343     int nstcalcenergy;
344     //! Group or verlet cutoffs
345     int cutoff_scheme;
346     //! Number of steps before pairlist is generated
347     int nstlist;
348     //! Number of steps after which center of mass motion is removed
349     int nstcomm;
350     //! Center of mass motion removal algorithm
351     int comm_mode;
352     //! Number of steps after which print to logfile
353     int nstlog;
354     //! Number of steps after which X is output
355     int nstxout;
356     //! Number of steps after which V is output
357     int nstvout;
358     //! Number of steps after which F is output
359     int nstfout;
360     //! Number of steps after which energies printed
361     int nstenergy;
362     //! Number of steps after which compressed trj (.xtc,.tng) is output
363     int nstxout_compressed;
364     //! Initial time (ps)
365     double init_t;
366     //! Time step (ps)
367     double delta_t;
368     //! Precision of x in compressed trajectory file
369     real x_compression_precision;
370     //! Requested fourier_spacing, when nk? not set
371     real fourier_spacing;
372     //! Number of k vectors in x dimension for fourier methods for long range electrost.
373     int nkx;
374     //! Number of k vectors in y dimension for fourier methods for long range electrost.
375     int nky;
376     //! Number of k vectors in z dimension for fourier methods for long range electrost.
377     int nkz;
378     //! Interpolation order for PME
379     int pme_order;
380     //! Real space tolerance for Ewald, determines the real/reciprocal space relative weight
381     real ewald_rtol;
382     //! Real space tolerance for LJ-Ewald
383     real ewald_rtol_lj;
384     //! Normal/3D ewald, or pseudo-2D LR corrections
385     int ewald_geometry;
386     //! Epsilon for PME dipole correction
387     real epsilon_surface;
388     //! Type of combination rule in LJ-PME
389     int ljpme_combination_rule;
390     //! Type of periodic boundary conditions
391     int ePBC;
392     //! Periodic molecules
393     bool bPeriodicMols;
394     //! Continuation run: starting state is correct (ie. constrained)
395     gmx_bool bContinuation;
396     //! Temperature coupling
397     int etc;
398     //! Interval in steps for temperature coupling
399     int nsttcouple;
400     //! Whether to print nose-hoover chains
401     gmx_bool bPrintNHChains;
402     //! Pressure coupling
403     int epc;
404     //! Pressure coupling type
405     int epct;
406     //! Interval in steps for pressure coupling
407     int nstpcouple;
408     //! Pressure coupling time (ps)
409     real tau_p;
410     //! Reference pressure (kJ/(mol nm^3))
411     tensor ref_p;
412     //! Compressibility ((mol nm^3)/kJ)
413     tensor compress;
414     //! How to scale absolute reference coordinates
415     int refcoord_scaling;
416     //! The COM of the posres atoms
417     rvec posres_com;
418     //! The B-state COM of the posres atoms
419     rvec posres_comB;
420     //! Random seed for Andersen thermostat (obsolete)
421     int andersen_seed;
422     //! Per atom pair energy drift tolerance (kJ/mol/ps/atom) for list buffer
423     real verletbuf_tol;
424     //! Short range pairlist cut-off (nm)
425     real rlist;
426     //! Radius for test particle insertion
427     real rtpi;
428     //! Type of electrostatics treatment
429     int coulombtype;
430     //! Modify the Coulomb interaction
431     int coulomb_modifier;
432     //! Coulomb switch range start (nm)
433     real rcoulomb_switch;
434     //! Coulomb cutoff (nm)
435     real rcoulomb;
436     //! Relative dielectric constant
437     real epsilon_r;
438     //! Relative dielectric constant of the RF
439     real epsilon_rf;
440     //! Always false (no longer supported)
441     bool implicit_solvent;
442     //! Type of Van der Waals treatment
443     int vdwtype;
444     //! Modify the Van der Waals interaction
445     int vdw_modifier;
446     //! Van der Waals switch range start (nm)
447     real rvdw_switch;
448     //! Van der Waals cutoff (nm)
449     real rvdw;
450     //! Perform Long range dispersion corrections
451     int eDispCorr;
452     //! Extension of the table beyond the cut-off, as well as the table length for 1-4 interac.
453     real tabext;
454     //! Tolerance for shake
455     real shake_tol;
456     //! Free energy calculations
457     int efep;
458     //! Data for the FEP state
459     t_lambda* fepvals;
460     //! Whether to do simulated tempering
461     gmx_bool bSimTemp;
462     //! Variables for simulated tempering
463     t_simtemp* simtempvals;
464     //! Whether expanded ensembles are used
465     gmx_bool bExpanded;
466     //! Expanded ensemble parameters
467     t_expanded* expandedvals;
468     //! Type of distance restraining
469     int eDisre;
470     //! Force constant for time averaged distance restraints
471     real dr_fc;
472     //! Type of weighting of pairs in one restraints
473     int eDisreWeighting;
474     //! Use combination of time averaged and instantaneous violations
475     gmx_bool bDisreMixed;
476     //! Frequency of writing pair distances to enx
477     int nstdisreout;
478     //! Time constant for memory function in disres
479     real dr_tau;
480     //! Force constant for orientational restraints
481     real orires_fc;
482     //! Time constant for memory function in orires
483     real orires_tau;
484     //! Frequency of writing tr(SD) to energy output
485     int nstorireout;
486     //! The stepsize for updating
487     real em_stepsize;
488     //! The tolerance
489     real em_tol;
490     //! Number of iterations for convergence of steepest descent in relax_shells
491     int niter;
492     //! Stepsize for directional minimization in relax_shells
493     real fc_stepsize;
494     //! Number of steps after which a steepest descents step is done while doing cg
495     int nstcgsteep;
496     //! Number of corrections to the Hessian to keep
497     int nbfgscorr;
498     //! Type of constraint algorithm
499     int eConstrAlg;
500     //! Order of the LINCS Projection Algorithm
501     int nProjOrder;
502     //! Warn if any bond rotates more than this many degrees
503     real LincsWarnAngle;
504     //! Number of iterations in the final LINCS step
505     int nLincsIter;
506     //! Use successive overrelaxation for shake
507     gmx_bool bShakeSOR;
508     //! Friction coefficient for BD (amu/ps)
509     real bd_fric;
510     //! Random seed for SD and BD
511     int64_t ld_seed;
512     //! The number of walls
513     int nwall;
514     //! The type of walls
515     int wall_type;
516     //! The potentail is linear for r<=wall_r_linpot
517     real wall_r_linpot;
518     //! The atom type for walls
519     int wall_atomtype[2];
520     //! Number density for walls
521     real wall_density[2];
522     //! Scaling factor for the box for Ewald
523     real wall_ewald_zfac;
524
525     /* COM pulling data */
526     //! Do we do COM pulling?
527     gmx_bool bPull;
528     //! The data for center of mass pulling
529     pull_params_t* pull;
530
531     /* AWH bias data */
532     //! Whether to use AWH biasing for PMF calculations
533     gmx_bool bDoAwh;
534     //! AWH biasing parameters
535     gmx::AwhParams* awhParams;
536
537     /* Enforced rotation data */
538     //! Whether to calculate enforced rotation potential(s)
539     gmx_bool bRot;
540     //! The data for enforced rotation potentials
541     t_rot* rot;
542
543     //! Whether to do ion/water position exchanges (CompEL)
544     int eSwapCoords;
545     //! Swap data structure.
546     t_swapcoords* swap;
547
548     //! Whether the tpr makes an interactive MD session possible.
549     gmx_bool bIMD;
550     //! Interactive molecular dynamics
551     t_IMD* imd;
552
553     //! Acceleration for viscosity calculation
554     real cos_accel;
555     //! Triclinic deformation velocities (nm/ps)
556     tensor deform;
557     /*! \brief User determined parameters */
558     /**@{*/
559     int  userint1;
560     int  userint2;
561     int  userint3;
562     int  userint4;
563     real userreal1;
564     real userreal2;
565     real userreal3;
566     real userreal4;
567     /**@}*/
568     //! Group options
569     t_grpopts opts;
570     //! QM/MM calculation
571     gmx_bool bQMMM;
572     //! Constraints on QM bonds
573     int QMconstraints;
574     //! Scheme: ONIOM or normal
575     int QMMMscheme;
576     //! Factor for scaling the MM charges in QM calc.
577     real scalefactor;
578
579     /* Fields for removed features go here (better caching) */
580     //! Whether AdResS is enabled - always false if a valid .tpr was read
581     gmx_bool bAdress;
582     //! Whether twin-range scheme is active - always false if a valid .tpr was read
583     gmx_bool useTwinRange;
584
585     //! KVT object that contains input parameters converted to the new style.
586     gmx::KeyValueTreeObject* params;
587
588     //! KVT for storing simulation parameters that are not part of the mdp file.
589     std::unique_ptr<gmx::KeyValueTreeObject> internalParameters;
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, gmx_bool bMDPformat);
637
638 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol);
639
640 void comp_pull_AB(FILE* fp, pull_params_t* pull, real ftol, real abstol);
641
642
643 gmx_bool inputrecDeform(const t_inputrec* ir);
644
645 gmx_bool inputrecDynamicBox(const t_inputrec* ir);
646
647 gmx_bool inputrecPreserveShape(const t_inputrec* ir);
648
649 gmx_bool inputrecNeedMutot(const t_inputrec* ir);
650
651 gmx_bool inputrecTwinRange(const t_inputrec* ir);
652
653 gmx_bool inputrecExclForces(const t_inputrec* ir);
654
655 gmx_bool inputrecNptTrotter(const t_inputrec* ir);
656
657 gmx_bool inputrecNvtTrotter(const t_inputrec* ir);
658
659 gmx_bool inputrecNphTrotter(const t_inputrec* ir);
660
661 /*! \brief Return true if the simulation is 2D periodic with two walls. */
662 bool inputrecPbcXY2Walls(const t_inputrec* ir);
663
664 /*! \brief Returns true for MD integator with T and/or P-coupling that supports
665  * calculating the conserved energy quantity.
666  */
667 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir);
668
669 /*! \brief Returns true when temperature is coupled or constant. */
670 bool integratorHasReferenceTemperature(const t_inputrec* ir);
671
672 /*! \brief Return the number of bounded dimensions
673  *
674  * \param[in] ir The input record with MD parameters
675  * \return the number of dimensions in which
676  * the coordinates of the particles are bounded, starting at X.
677  */
678 int inputrec2nboundeddim(const t_inputrec* ir);
679
680 /*! \brief Returns the number of degrees of freedom in center of mass motion
681  *
682  * \param[in] ir  The inputrec structure
683  * \return the number of degrees of freedom of the center of mass
684  */
685 int ndof_com(const t_inputrec* ir);
686
687 /*! \brief Returns the maximum reference temperature over all coupled groups
688  *
689  * Returns 0 for energy minimization and normal mode computation.
690  * Returns -1 for MD without temperature coupling.
691  *
692  * \param[in] ir  The inputrec structure
693  */
694 real maxReferenceTemperature(const t_inputrec& ir);
695
696 /*! \brief Returns whether there is an Ewald surface contribution
697  */
698 bool haveEwaldSurfaceContribution(const t_inputrec& ir);
699
700 #endif /* GMX_MDTYPES_INPUTREC_H */