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