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