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