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