95cfc402cd00bb5cc3177bea6bfa1dc5c0218e38
[alexxy/gromacs.git] / src / gromacs / legacyheaders / types / 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, 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 _inputrec_h_
38 #define _inputrec_h_
39
40 #include <stdio.h>
41
42 #include "gromacs/legacyheaders/types/enums.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/swap/enums.h"
45
46 #ifdef __cplusplus
47 extern "C" {
48 #endif
49
50
51 typedef struct {
52     int   n;    /* Number of terms                              */
53     real *a;    /* Coeffients (V / nm )                     */
54     real *phi;  /* Phase angles                                 */
55 } t_cosines;
56
57 typedef struct {
58     real E0;            /* Field strength (V/nm)                        */
59     real omega;         /* Frequency (1/ps)                             */
60     real t0;            /* Centre of the Gaussian pulse (ps)            */
61     real sigma;         /* Width of the Gaussian pulse (FWHM) (ps)      */
62 } t_efield;
63
64 #define EGP_EXCL  (1<<0)
65 #define EGP_TABLE (1<<1)
66
67 typedef struct {
68     int       ngtc;           /* # T-Coupl groups                        */
69     int       nhchainlength;  /* # of nose-hoover chains per group       */
70     int       ngacc;          /* # Accelerate groups                     */
71     int       ngfrz;          /* # Freeze groups                         */
72     int       ngener;         /* # Ener groups                      */
73     real     *nrdf;           /* Nr of degrees of freedom in a group        */
74     real     *ref_t;          /* Coupling temperature   per group   */
75     int      *annealing;      /* No/simple/periodic SA for each group    */
76     int      *anneal_npoints; /* Number of annealing time points per grp */
77     real    **anneal_time;    /* For ea. group: Time points              */
78     real    **anneal_temp;    /* For ea. grp: Temperature at these times */
79                               /* Final temp after all intervals is ref_t */
80     real     *tau_t;          /* Tau coupling time              */
81     rvec     *acc;            /* Acceleration per group             */
82     ivec     *nFreeze;        /* Freeze the group in each direction ?    */
83     int      *egp_flags;      /* Exclusions/tables of energy group pairs */
84
85     /* QMMM stuff */
86     int          ngQM;         /* nr of QM groups                              */
87     int         *QMmethod;     /* Level of theory in the QM calculation        */
88     int         *QMbasis;      /* Basisset in the QM calculation               */
89     int         *QMcharge;     /* Total charge in the QM region                */
90     int         *QMmult;       /* Spin multiplicicty in the QM region          */
91     gmx_bool    *bSH;          /* surface hopping (diabatic hop only)          */
92     int         *CASorbitals;  /* number of orbiatls in the active space       */
93     int         *CASelectrons; /* number of electrons in the active space      */
94     real        *SAon;         /* at which gap (A.U.) the SA is switched on    */
95     real        *SAoff;
96     int         *SAsteps;      /* in how many steps SA goes from 1-1 to 0.5-0.5*/
97     gmx_bool    *bOPT;
98     gmx_bool    *bTS;
99 } t_grpopts;
100
101 enum {
102     epgrppbcNONE, epgrppbcREFAT, epgrppbcCOS
103 };
104
105 typedef struct {
106     int         nat;        /* Number of atoms in the pull group */
107     atom_id    *ind;        /* The global atoms numbers */
108     int         nat_loc;    /* Number of local pull atoms */
109     int         nalloc_loc; /* Allocation size for ind_loc and weight_loc */
110     atom_id    *ind_loc;    /* Local pull indices */
111     int         nweight;    /* The number of weights (0 or nat) */
112     real       *weight;     /* Weights (use all 1 when weight==NULL) */
113     real       *weight_loc; /* Weights for the local indices */
114     int         epgrppbc;   /* The type of pbc for this pull group, see enum above */
115     atom_id     pbcatom;    /* The reference atom for pbc (global number) */
116
117     /* Variables not present in mdp, but used at run time */
118     gmx_bool    bCalcCOM;   /* Calculate COM? Not if only used as cylinder group */
119     real        mwscale;    /* mass*weight scaling factor 1/sum w m */
120     real        wscale;     /* scaling factor for the weights: sum w m/sum w w m */
121     real        invtm;      /* inverse total mass of the group: 1/wscale sum w m */
122     dvec       *mdw;        /* mass*gradient(weight) for atoms */
123     double     *dv;         /* distance to the other group along vec */
124     dvec        x;          /* center of mass before update */
125     dvec        xp;         /* center of mass after update before constraining */
126 } t_pull_group;
127
128 typedef struct {
129     int         group[4];   /* The pull groups, index in group in t_pull */
130     int         eType;      /* The pull type: umbrella, constraint, ... */
131     int         eGeom;      /* The pull geometry */
132     ivec        dim;        /* Used to select components for constraint */
133     rvec        origin;     /* The origin for the absolute reference */
134     rvec        vec;        /* The pull vector, direction or position */
135     gmx_bool    bStart;     /* Set init based on the initial structure */
136     real        init;       /* Initial reference displacement */
137     real        rate;       /* Rate of motion (nm/ps) */
138     real        k;          /* force constant */
139     real        kB;         /* force constant for state B */
140
141     /* Variables not present in mdp, but used at run time */
142     double      value_ref;  /* The reference value, usually init+rate*t */
143     double      value;      /* The current value of the coordinate */
144     dvec        dr;         /* The distance from the reference group */
145     double      vec_len;    /* Length of vec for direction-relative */
146     dvec        ffrad;      /* conversion factor from vec to radial force */
147     double      cyl_dev;    /* The deviation from the reference position */
148     double      f_scal;     /* Scalar force for directional pulling */
149     dvec        f;          /* force due to the pulling/constraining */
150 } t_pull_coord;
151
152 typedef struct {
153     int   eSimTempScale; /* simulated temperature scaling; linear or exponential */
154     real  simtemp_low;   /* the low temperature for simulated tempering  */
155     real  simtemp_high;  /* the high temperature for simulated tempering */
156     real *temperatures;  /* the range of temperatures used for simulated tempering */
157 } t_simtemp;
158
159 typedef struct {
160     int    nstdhdl;                 /* The frequency for calculating dhdl           */
161     double init_lambda;             /* fractional value of lambda (usually will use
162                                        init_fep_state, this will only be for slow growth,
163                                        and for legacy free energy code. Only has a
164                                        valid value if positive)   */
165     int      init_fep_state;        /* the initial number of the state                   */
166     double   delta_lambda;          /* change of lambda per time step (fraction of (0.1) */
167     int      edHdLPrintEnergy;      /* print no, total or potential energies in dhdl    */
168     int      n_lambda;              /* The number of foreign lambda points               */
169     double **all_lambda;            /* The array of all lambda values                    */
170     int      lambda_neighbors;      /* The number of neighboring lambda states to
171                                        calculate the energy for in up and down directions
172                                        (-1 for all) */
173     int      lambda_start_n;        /* The first lambda to calculate energies for */
174     int      lambda_stop_n;         /* The last lambda +1 to calculate energies for */
175     real     sc_alpha;              /* free energy soft-core parameter                   */
176     int      sc_power;              /* lambda power for soft-core interactions           */
177     real     sc_r_power;            /* r power for soft-core interactions                */
178     real     sc_sigma;              /* free energy soft-core sigma when c6 or c12=0      */
179     real     sc_sigma_min;          /* free energy soft-core sigma for ?????             */
180     gmx_bool bScCoul;               /* use softcore for the coulomb portion as well (default FALSE) */
181     gmx_bool separate_dvdl[efptNR]; /* whether to print the dvdl term associated with
182                                        this term; if it is not specified as separate,
183                                        it is lumped with the FEP term */
184     int    separate_dhdl_file;      /* whether to write a separate dhdl.xvg file
185                                        note: NOT a gmx_bool, but an enum */
186     int    dhdl_derivatives;        /* whether to calculate+write dhdl derivatives
187                                        note: NOT a gmx_bool, but an enum */
188     int    dh_hist_size;            /* The maximum table size for the dH histogram */
189     double dh_hist_spacing;         /* The spacing for the dH histogram */
190 } t_lambda;
191
192 typedef struct {
193     int      nstexpanded;         /* The frequency of expanded ensemble state changes */
194     int      elamstats;           /* which type of move updating do we use for lambda monte carlo (or no for none) */
195     int      elmcmove;            /* what move set will be we using for state space moves */
196     int      elmceq;              /* the method we use to decide of we have equilibrated the weights */
197     int      equil_n_at_lam;      /* the minumum number of samples at each lambda for deciding whether we have reached a minimum */
198     real     equil_wl_delta;      /* WL delta at which we stop equilibrating weights */
199     real     equil_ratio;         /* use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating */
200     int      equil_steps;         /* after equil_steps steps we stop equilibrating the weights */
201     int      equil_samples;       /* after equil_samples total samples (steps/nstfep), we stop equilibrating the weights */
202     int      lmc_seed;            /* random number seed for lambda mc switches */
203     gmx_bool minvar;              /* whether to use minumum variance weighting */
204     int      minvarmin;           /* the number of samples needed before kicking into minvar routine */
205     real     minvar_const;        /* the offset for the variance in MinVar */
206     int      c_range;             /* range of cvalues used for BAR */
207     gmx_bool bSymmetrizedTMatrix; /* whether to print symmetrized matrices */
208     int      nstTij;              /* How frequently to print the transition matrices */
209     int      lmc_repeats;         /* number of repetitions in the MC lambda jumps */  /*MRS -- VERIFY THIS */
210     int      lmc_forced_nstart;   /* minimum number of samples for each state before free sampling */ /* MRS -- VERIFY THIS! */
211     int      gibbsdeltalam;       /* distance in lambda space for the gibbs interval */
212     real     wl_scale;            /* scaling factor for wang-landau */
213     real     wl_ratio;            /* ratio between largest and smallest number for freezing the weights */
214     real     init_wl_delta;       /* starting delta for wang-landau */
215     gmx_bool bWLoneovert;         /* use one over t convergece for wang-landau when the delta get sufficiently small */
216     gmx_bool bInit_weights;       /* did we initialize the weights? TODO: REMOVE FOR 5.0, no longer needed with new logic */
217     real     mc_temp;             /* To override the main temperature, or define it if it's not defined */
218     real    *init_lambda_weights; /* user-specified initial weights to start with  */
219 } t_expanded;
220
221 typedef struct {
222     int            ngroup;         /* number of pull groups */
223     int            ncoord;         /* number of pull coordinates */
224     real           cylinder_r;     /* radius of cylinder for dynamic COM */
225     real           constr_tol;     /* absolute tolerance for constraints in (nm) */
226     gmx_bool       bPrintCOM1;     /* Print coordinates of COM 1 for each coord */
227     gmx_bool       bPrintCOM2;     /* Print coordinates of COM 2 for each coord */
228     gmx_bool       bPrintRefValue; /* Print the reference value for each coord */
229     gmx_bool       bPrintComp;     /* Print cartesian components for each coord with geometry=distance */
230     int            nstxout;        /* Output frequency for pull x */
231     int            nstfout;        /* Output frequency for pull f */
232     int            ePBC;           /* the boundary conditions */
233     int            npbcdim;        /* do pbc in dims 0 <= dim < npbcdim */
234     gmx_bool       bRefAt;         /* do we need reference atoms for a group COM ? */
235     int            cosdim;         /* dimension for cosine weighting, -1 if none */
236     t_pull_group  *group;          /* groups to pull/restrain/etc/ */
237     t_pull_coord  *coord;          /* the pull coordinates */
238
239     /* Variables not present in mdp, but used at run time */
240     gmx_bool       bPotential;   /* Are there coordinates with potential? */
241     gmx_bool       bConstraint;  /* Are there constrained coordinates? */
242     gmx_bool       bCylinder;    /* Is group 0 a cylinder group? */
243     t_pull_group  *dyna;         /* dynamic groups for use with local constraints */
244     gmx_bool       bSetPBCatoms; /* Do we need to set x_pbc for the groups? */
245
246     rvec          *rbuf;         /* COM calculation buffer */
247     dvec          *dbuf;         /* COM calculation buffer */
248     double        *dbuf_cyl;     /* cylinder ref. groups COM calculation buffer */
249
250     FILE          *out_x;        /* output file for pull data */
251     FILE          *out_f;        /* output file for pull data */
252 } t_pull;
253
254
255 /* Abstract types for enforced rotation only defined in pull_rotation.c       */
256 typedef struct gmx_enfrot *gmx_enfrot_t;
257 typedef struct gmx_enfrotgrp *gmx_enfrotgrp_t;
258
259 typedef struct {
260     int         eType;             /* Rotation type for this group                  */
261     int         bMassW;            /* Use mass-weighed positions?                   */
262     int         nat;               /* Number of atoms in the group                  */
263     atom_id    *ind;               /* The global atoms numbers                      */
264     rvec       *x_ref;             /* The reference positions                       */
265     rvec        vec;               /* The normalized rotation vector                */
266     real        rate;              /* Rate of rotation (degree/ps)                  */
267     real        k;                 /* Force constant (kJ/(mol nm^2)                 */
268     rvec        pivot;             /* Pivot point of rotation axis (nm)             */
269     int         eFittype;          /* Type of fit to determine actual group angle   */
270     int         PotAngle_nstep;    /* Number of angles around the reference angle
271                                       for which the rotation potential is also
272                                       evaluated (for fit type 'potential' only)     */
273     real            PotAngle_step; /* Distance between two angles in degrees (for
274                                       fit type 'potential' only)                    */
275     real            slab_dist;     /* Slab distance (nm)                            */
276     real            min_gaussian;  /* Minimum value the gaussian must have so that
277                                       the force is actually evaluated               */
278     real            eps;           /* Additive constant for radial motion2 and
279                                       flexible2 potentials (nm^2)                   */
280     gmx_enfrotgrp_t enfrotgrp;     /* Stores non-inputrec rotation data per group   */
281 } t_rotgrp;
282
283 typedef struct {
284     int          ngrp;       /* Number of rotation groups                     */
285     int          nstrout;    /* Output frequency for main rotation outfile    */
286     int          nstsout;    /* Output frequency for per-slab data            */
287     t_rotgrp    *grp;        /* Groups to rotate                              */
288     gmx_enfrot_t enfrot;     /* Stores non-inputrec enforced rotation data    */
289 } t_rot;
290
291 /* Abstract type for IMD only defined in IMD.c */
292 typedef struct gmx_IMD *t_gmx_IMD;
293
294 typedef struct {
295     int         nat;         /* Number of interactive atoms                   */
296     atom_id    *ind;         /* The global indices of the interactive atoms   */
297     t_gmx_IMD   setup;       /* Stores non-inputrec IMD data                  */
298 } t_IMD;
299
300 /* Abstract types for position swapping only defined in swapcoords.c */
301 typedef struct t_swap *gmx_swapcoords_t;
302
303 typedef struct {
304     int              nstswap;           /* Every how many steps a swap is attempted?    */
305     int              nat;               /* Number of atoms in the ion group             */
306     int              nat_split[2];      /* Number of atoms in the split group           */
307     int              nat_sol;           /* Number of atoms in the solvent group         */
308     atom_id         *ind;               /* The global ion group atoms numbers           */
309     atom_id         *ind_split[2];      /* Split groups for compartment partitioning    */
310     atom_id         *ind_sol;           /* The global solvent group atom numbers        */
311     gmx_bool         massw_split[2];    /* Use mass-weighted positions in split group?  */
312     real             cyl0r, cyl1r;      /* Split cylinders defined by radius, upper and */
313     real             cyl0u, cyl1u;      /* ... lower extension. The split cylinders de- */
314     real             cyl0l, cyl1l;      /* ... fine the channels and are each anchored  */
315                                         /* ... in the center of the split group         */
316     int              nanions[eCompNR];  /* Requested number of anions and               */
317     int              nAverage;          /* Coupling constant (nr of swap attempt steps) */
318     real             threshold;         /* Ion counts may deviate from the requested
319                                            values by +-threshold before a swap is done  */
320     int              ncations[eCompNR]; /* ... cations for both compartments            */
321     gmx_swapcoords_t si_priv;           /* swap private data accessible in
322                                          * swapcoords.c                                 */
323 } t_swapcoords;
324
325
326 typedef struct {
327     int      type;           /* type of AdResS simulation                    */
328     gmx_bool bnew_wf;        /* enable new AdResS weighting function         */
329     gmx_bool bchempot_dx;    /*true:interaction table format input is F=-dmu/dx   false: dmu_dwp  */
330     gmx_bool btf_full_box;   /* true: appy therm force everywhere in the box according to table false: only in hybrid region */
331     real     const_wf;       /* value of weighting function for eAdressConst */
332     real     ex_width;       /* center of the explicit zone                  */
333     real     hy_width;       /* width of the hybrid zone                     */
334     int      icor;           /* type of interface correction                 */
335     int      site;           /* AdResS CG site location                      */
336     rvec     refs;           /* Coordinates for AdResS reference             */
337     real     ex_forcecap;    /* in the hybrid zone, cap forces large then this to adress_ex_forcecap */
338     gmx_bool do_hybridpairs; /* If true pair interaction forces are also scaled in an adress way*/
339
340     int    * tf_table_index; /* contains mapping of energy group index -> i-th adress tf table*/
341     int      n_tf_grps;
342     int     *group_explicit;
343     int      n_energy_grps;
344 } t_adress;
345
346 typedef struct {
347     int             eI;                      /* Integration method                 */
348     gmx_int64_t     nsteps;                  /* number of steps to be taken                     */
349     int             simulation_part;         /* Used in checkpointing to separate chunks */
350     gmx_int64_t     init_step;               /* start at a stepcount >0 (used w. convert-tpr)    */
351     int             nstcalcenergy;           /* frequency of energy calc. and T/P coupl. upd.   */
352     int             cutoff_scheme;           /* group or verlet cutoffs     */
353     int             ns_type;                 /* which ns method should we use?               */
354     int             nstlist;                 /* number of steps before pairlist is generated    */
355     int             ndelta;                  /* number of cells per rlong                       */
356     int             nstcomm;                 /* number of steps after which center of mass      */
357                                              /* motion is removed                               */
358     int             comm_mode;               /* Center of mass motion removal algorithm      */
359     int             nstlog;                  /* number of steps after which print to logfile    */
360     int             nstxout;                 /* number of steps after which X is output */
361     int             nstvout;                 /* id. for V                                       */
362     int             nstfout;                 /* id. for F                                       */
363     int             nstenergy;               /* number of steps after which energies printed */
364     int             nstxout_compressed;      /* id. for compressed trj (.xtc,.tng)           */
365     double          init_t;                  /* initial time (ps)              */
366     double          delta_t;                 /* time step (ps)                          */
367     real            x_compression_precision; /* precision of x in compressed trajectory file */
368     real            fourier_spacing;         /* requested fourier_spacing, when nk? not set  */
369     int             nkx, nky, nkz;           /* number of k vectors in each spatial dimension*/
370                                              /* for fourier methods for long range electrost.*/
371     int             pme_order;               /* interpolation order for PME                  */
372     real            ewald_rtol;              /* Real space tolerance for Ewald, determines   */
373                                              /* the real/reciprocal space relative weight    */
374     real            ewald_rtol_lj;           /* Real space tolerance for LJ-Ewald            */
375     int             ewald_geometry;          /* normal/3d ewald, or pseudo-2d LR corrections */
376     real            epsilon_surface;         /* Epsilon for PME dipole correction            */
377     int             ljpme_combination_rule;  /* Type of combination rule in LJ-PME          */
378     int             ePBC;                    /* Type of periodic boundary conditions            */
379     int             bPeriodicMols;           /* Periodic molecules                           */
380     gmx_bool        bContinuation;           /* Continuation run: starting state is correct     */
381     int             etc;                     /* temperature coupling               */
382     int             nsttcouple;              /* interval in steps for temperature coupling   */
383     gmx_bool        bPrintNHChains;          /* whether to print nose-hoover chains        */
384     int             epc;                     /* pressure coupling                            */
385     int             epct;                    /* pressure coupling type                  */
386     int             nstpcouple;              /* interval in steps for pressure coupling      */
387     real            tau_p;                   /* pressure coupling time (ps)                     */
388     tensor          ref_p;                   /* reference pressure (kJ/(mol nm^3))              */
389     tensor          compress;                /* compressability ((mol nm^3)/kJ)        */
390     int             refcoord_scaling;        /* How to scale absolute reference coordinates  */
391     rvec            posres_com;              /* The COM of the posres atoms                  */
392     rvec            posres_comB;             /* The B-state COM of the posres atoms          */
393     int             andersen_seed;           /* Random seed for Andersen thermostat (obsolete) */
394     real            verletbuf_tol;           /* Per atom pair energy drift tolerance (kJ/mol/ps/atom) for list buffer  */
395     real            rlist;                   /* short range pairlist cut-off (nm)               */
396     real            rlistlong;               /* long range pairlist cut-off (nm)                */
397     int             nstcalclr;               /* Frequency of evaluating direct space long-range interactions */
398     real            rtpi;                    /* Radius for test particle insertion           */
399     int             coulombtype;             /* Type of electrostatics treatment             */
400     int             coulomb_modifier;        /* Modify the Coulomb interaction              */
401     real            rcoulomb_switch;         /* Coulomb switch range start (nm)         */
402     real            rcoulomb;                /* Coulomb cutoff (nm)                             */
403     real            epsilon_r;               /* relative dielectric constant                 */
404     real            epsilon_rf;              /* relative dielectric constant of the RF       */
405     int             implicit_solvent;        /* No (=explicit water), or GBSA solvent models */
406     int             gb_algorithm;            /* Algorithm to use for calculation Born radii  */
407     int             nstgbradii;              /* Frequency of updating Generalized Born radii */
408     real            rgbradii;                /* Cutoff for GB radii calculation              */
409     real            gb_saltconc;             /* Salt concentration (M) for GBSA models       */
410     real            gb_epsilon_solvent;      /* dielectric coeff. of implicit solvent     */
411     real            gb_obc_alpha;            /* 1st scaling factor for Bashford-Case GB      */
412     real            gb_obc_beta;             /* 2nd scaling factor for Bashford-Case GB      */
413     real            gb_obc_gamma;            /* 3rd scaling factor for Bashford-Case GB      */
414     real            gb_dielectric_offset;    /* Dielectric offset for Still/HCT/OBC     */
415     int             sa_algorithm;            /* Algorithm for SA part of GBSA                */
416     real            sa_surface_tension;      /* Energy factor for SA part of GBSA */
417     int             vdwtype;                 /* Type of Van der Waals treatment              */
418     int             vdw_modifier;            /* Modify the VdW interaction                   */
419     real            rvdw_switch;             /* Van der Waals switch range start (nm)        */
420     real            rvdw;                    /* Van der Waals cutoff (nm)               */
421     int             eDispCorr;               /* Perform Long range dispersion corrections    */
422     real            tabext;                  /* Extension of the table beyond the cut-off,   *
423                                               * as well as the table length for 1-4 interac. */
424     real            shake_tol;               /* tolerance for shake                             */
425     int             efep;                    /* free energy calculations                     */
426     t_lambda       *fepvals;                 /* Data for the FEP state                       */
427     gmx_bool        bSimTemp;                /* Whether to do simulated tempering            */
428     t_simtemp      *simtempvals;             /* Variables for simulated tempering            */
429     gmx_bool        bExpanded;               /* Whether expanded ensembles are used          */
430     t_expanded     *expandedvals;            /* Expanded ensemble parameters              */
431     int             eDisre;                  /* Type of distance restraining                 */
432     real            dr_fc;                   /* force constant for ta_disre                     */
433     int             eDisreWeighting;         /* type of weighting of pairs in one restraints    */
434     gmx_bool        bDisreMixed;             /* Use comb of time averaged and instan. viol's    */
435     int             nstdisreout;             /* frequency of writing pair distances to enx   */
436     real            dr_tau;                  /* time constant for memory function in disres    */
437     real            orires_fc;               /* force constant for orientational restraints  */
438     real            orires_tau;              /* time constant for memory function in orires    */
439     int             nstorireout;             /* frequency of writing tr(SD) to enx           */
440     real            em_stepsize;             /* The stepsize for updating                       */
441     real            em_tol;                  /* The tolerance                           */
442     int             niter;                   /* Number of iterations for convergence of      */
443                                              /* steepest descent in relax_shells             */
444     real            fc_stepsize;             /* Stepsize for directional minimization        */
445                                              /* in relax_shells                              */
446     int             nstcgsteep;              /* number of steps after which a steepest       */
447                                              /* descents step is done while doing cg         */
448     int             nbfgscorr;               /* Number of corrections to the hessian to keep */
449     int             eConstrAlg;              /* Type of constraint algorithm                 */
450     int             nProjOrder;              /* Order of the LINCS Projection Algorithm      */
451     real            LincsWarnAngle;          /* If bond rotates more than %g degrees, warn   */
452     int             nLincsIter;              /* Number of iterations in the final Lincs step */
453     gmx_bool        bShakeSOR;               /* Use successive overrelaxation for shake      */
454     real            bd_fric;                 /* Friction coefficient for BD (amu/ps)         */
455     gmx_int64_t     ld_seed;                 /* Random seed for SD and BD                    */
456     int             nwall;                   /* The number of walls                          */
457     int             wall_type;               /* The type of walls                            */
458     real            wall_r_linpot;           /* The potentail is linear for r<=wall_r_linpot */
459     int             wall_atomtype[2];        /* The atom type for walls                      */
460     real            wall_density[2];         /* Number density for walls                     */
461     real            wall_ewald_zfac;         /* Scaling factor for the box for Ewald         */
462     gmx_bool        bPull;                   /* Do we do COM pulling?                        */
463     t_pull         *pull;                    /* The data for center of mass pulling          */
464     gmx_bool        bRot;                    /* Calculate enforced rotation potential(s)?    */
465     t_rot          *rot;                     /* The data for enforced rotation potentials    */
466     int             eSwapCoords;             /* Do ion/water position exchanges (CompEL)?    */
467     t_swapcoords   *swap;
468     gmx_bool        bIMD;                    /* Allow interactive MD sessions for this .tpr? */
469     t_IMD          *imd;                     /* Interactive molecular dynamics               */
470     real            cos_accel;               /* Acceleration for viscosity calculation       */
471     tensor          deform;                  /* Triclinic deformation velocities (nm/ps)     */
472     int             userint1;                /* User determined parameters                   */
473     int             userint2;
474     int             userint3;
475     int             userint4;
476     real            userreal1;
477     real            userreal2;
478     real            userreal3;
479     real            userreal4;
480     t_grpopts       opts;          /* Group options                             */
481     t_cosines       ex[DIM];       /* Electric field stuff      (spatial part)          */
482     t_cosines       et[DIM];       /* Electric field stuff      (time part)             */
483     gmx_bool        bQMMM;         /* QM/MM calculation                            */
484     int             QMconstraints; /* constraints on QM bonds                      */
485     int             QMMMscheme;    /* Scheme: ONIOM or normal                      */
486     real            scalefactor;   /* factor for scaling the MM charges in QM calc.*/
487                                    /* parameter needed for AdResS simulation       */
488     gmx_bool        bAdress;       /* Is AdResS enabled ? */
489     t_adress       *adress;        /* The data for adress simulations */
490 } t_inputrec;
491
492 #define DEFORM(ir) ((ir).deform[XX][XX] != 0 || (ir).deform[YY][YY] != 0 || (ir).deform[ZZ][ZZ] != 0 || (ir).deform[YY][XX] != 0 || (ir).deform[ZZ][XX] != 0 || (ir).deform[ZZ][YY] != 0)
493
494 #define DYNAMIC_BOX(ir) ((ir).epc != epcNO || (ir).eI == eiTPI || DEFORM(ir))
495
496 #define PRESERVE_SHAPE(ir) ((ir).epc != epcNO && (ir).deform[XX][XX] == 0 && ((ir).epct == epctISOTROPIC || (ir).epct == epctSEMIISOTROPIC))
497
498 #define NEED_MUTOT(ir) (((ir).coulombtype == eelEWALD || EEL_PME((ir).coulombtype)) && ((ir).ewald_geometry == eewg3DC || (ir).epsilon_surface != 0))
499
500 #define IR_TWINRANGE(ir) ((ir).rlist > 0 && ((ir).rlistlong == 0 || (ir).rlistlong > (ir).rlist))
501
502 #define IR_ELEC_FIELD(ir) ((ir).ex[XX].n > 0 || (ir).ex[YY].n > 0 || (ir).ex[ZZ].n > 0)
503
504 #define IR_EXCL_FORCES(ir) (EEL_FULL((ir).coulombtype) || (EEL_RF((ir).coulombtype) && (ir).coulombtype != eelRF_NEC) || (ir).implicit_solvent != eisNO)
505 /* use pointer definitions of ir here, since that's what's usually used in the code */
506 #define IR_NPT_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && (((ir)->epc == epcMTTK) && ((ir)->etc == etcNOSEHOOVER)))
507
508 #define IR_NVT_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && ((!((ir)->epc == epcMTTK)) && ((ir)->etc == etcNOSEHOOVER)))
509
510 #define IR_NPH_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && (((ir)->epc == epcMTTK) && (!(((ir)->etc == etcNOSEHOOVER)))))
511
512 #ifdef __cplusplus
513 }
514 #endif
515
516
517 #endif