added mdp parameter for FE dH tables
[gromacs.git] / include / types / inputrec.h
blob869578dcbf815fb4ef7b1f77cd1eac3435a6c858
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GRoups of Organic Molecules in ACtion for Science
35 #ifndef _inputrec_h_
36 #define _inputrec_h_
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
44 #include "simple.h"
46 #ifdef __cplusplus
47 extern "C" {
48 #endif
51 typedef struct {
52 int n; /* Number of terms */
53 real *a; /* Coeffients (V / nm ) */
54 real *phi; /* Phase angles */
55 } t_cosines;
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;
64 #define EGP_EXCL (1<<0)
65 #define EGP_TABLE (1<<1)
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 */
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 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 bool *bOPT;
98 bool *bTS;
99 } t_grpopts;
101 enum { epgrppbcNONE, epgrppbcREFAT, epgrppbcCOS };
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 rvec vec; /* The pull vector, direction or position */
115 rvec init; /* Initial reference displacement */
116 real rate; /* Rate of motion (nm/ps) */
117 real k; /* force constant */
118 real kB; /* force constant for state B */
119 real wscale; /* scaling factor for the weights: sum w m/sum w w m */
120 real invtm; /* inverse total mass of the group: 1/wscale sum w m */
121 dvec x; /* center of mass before update */
122 dvec xp; /* center of mass after update before constraining */
123 dvec dr; /* The distance from the reference group */
124 double f_scal; /* Scalar force for directional pulling */
125 dvec f; /* force due to the pulling/constraining */
126 } t_pullgrp;
128 typedef struct {
129 int ngrp; /* number of groups */
130 int eGeom; /* pull geometry */
131 ivec dim; /* used to select components for constraint */
132 real cyl_r1; /* radius of cylinder for dynamic COM */
133 real cyl_r0; /* radius of cylinder including switch length */
134 real constr_tol; /* absolute tolerance for constraints in (nm) */
135 int nstxout; /* Output frequency for pull x */
136 int nstfout; /* Output frequency for pull f */
137 int ePBC; /* the boundary conditions */
138 int npbcdim; /* do pbc in dims 0 <= dim < npbcdim */
139 bool bRefAt; /* do we need reference atoms for a group COM ? */
140 int cosdim; /* dimension for cosine weighting, -1 if none */
141 bool bVirial; /* do we need to add the pull virial? */
142 t_pullgrp *grp; /* groups to pull/restrain/etc/ */
143 t_pullgrp *dyna; /* dynamic groups for use with local constraints */
144 rvec *rbuf; /* COM calculation buffer */
145 dvec *dbuf; /* COM calculation buffer */
146 double *dbuf_cyl; /* cylinder ref. groups COM calculation buffer */
148 FILE *out_x; /* output file for pull data */
149 FILE *out_f; /* output file for pull data */
150 } t_pull;
152 typedef struct {
153 int eI; /* Integration method */
154 gmx_large_int_t nsteps; /* number of steps to be taken */
155 int simulation_part; /* Used in checkpointing to separate chunks */
156 gmx_large_int_t init_step; /* start at a stepcount >0 (used w. tpbconv) */
157 int nstcalcenergy; /* fequency of energy calc. and T/P coupl. upd. */
158 int ns_type; /* which ns method should we use? */
159 int nstlist; /* number of steps before pairlist is generated */
160 int ndelta; /* number of cells per rlong */
161 int nstcomm; /* number of steps after which center of mass */
162 /* motion is removed */
163 int comm_mode; /* Center of mass motion removal algorithm */
164 int nstcheckpoint; /* checkpointing frequency */
165 int nstlog; /* number of steps after which print to logfile */
166 int nstxout; /* number of steps after which X is output */
167 int nstvout; /* id. for V */
168 int nstfout; /* id. for F */
169 int nstenergy; /* number of steps after which energies printed */
170 int nstxtcout; /* id. for compressed trj (.xtc) */
171 double init_t; /* initial time (ps) */
172 double delta_t; /* time step (ps) */
173 real xtcprec; /* precision of xtc file */
174 int nkx,nky,nkz; /* number of k vectors in each spatial dimension*/
175 /* for fourier methods for long range electrost.*/
176 int pme_order; /* interpolation order for PME */
177 real ewald_rtol; /* Real space tolerance for Ewald, determines */
178 /* the real/reciprocal space relative weight */
179 int ewald_geometry; /* normal/3d ewald, or pseudo-2d LR corrections */
180 real epsilon_surface; /* Epsilon for PME dipole correction */
181 bool bOptFFT; /* optimize the fft plan at start */
182 int ePBC; /* Type of periodic boundary conditions */
183 int bPeriodicMols; /* Periodic molecules */
184 bool bContinuation; /* Continuation run: starting state is correct */
185 int etc; /* temperature coupling */
186 int nsttcouple; /* interval in steps for temperature coupling */
187 int epc; /* pressure coupling */
188 int epct; /* pressure coupling type */
189 int nstpcouple; /* interval in steps for pressure coupling */
190 real tau_p; /* pressure coupling time (ps) */
191 tensor ref_p; /* reference pressure (kJ/(mol nm^3)) */
192 tensor compress; /* compressability ((mol nm^3)/kJ) */
193 int refcoord_scaling;/* How to scale absolute reference coordinates */
194 rvec posres_com; /* The COM of the posres atoms */
195 rvec posres_comB; /* The B-state COM of the posres atoms */
196 int andersen_seed; /* Random seed for Andersen thermostat. */
197 real rlist; /* short range pairlist cut-off (nm) */
198 real rlistlong; /* long range pairlist cut-off (nm) */
199 real rtpi; /* Radius for test particle insertion */
200 int coulombtype; /* Type of electrostatics treatment */
201 real rcoulomb_switch; /* Coulomb switch range start (nm) */
202 real rcoulomb; /* Coulomb cutoff (nm) */
203 real epsilon_r; /* relative dielectric constant */
204 real epsilon_rf; /* relative dielectric constant of the RF */
205 int implicit_solvent;/* No (=explicit water), or GBSA solvent models */
206 int gb_algorithm; /* Algorithm to use for calculation Born radii */
207 int nstgbradii; /* Frequency of updating Generalized Born radii */
208 real rgbradii; /* Cutoff for GB radii calculation */
209 real gb_saltconc; /* Salt concentration (M) for GBSA models */
210 real gb_epsilon_solvent; /* dielectric coeff. of implicit solvent */
211 real gb_obc_alpha; /* 1st scaling factor for Bashford-Case GB */
212 real gb_obc_beta; /* 2nd scaling factor for Bashford-Case GB */
213 real gb_obc_gamma; /* 3rd scaling factor for Bashford-Case GB */
214 real gb_dielectric_offset; /* Dielectric offset for Still/HCT/OBC */
215 int sa_algorithm; /* Algorithm for SA part of GBSA */
216 real sa_surface_tension; /* Energy factor for SA part of GBSA */
217 int vdwtype; /* Type of Van der Waals treatment */
218 real rvdw_switch; /* Van der Waals switch range start (nm) */
219 real rvdw; /* Van der Waals cutoff (nm) */
220 int eDispCorr; /* Perform Long range dispersion corrections */
221 real tabext; /* Extension of the table beyond the cut-off, *
222 * as well as the table length for 1-4 interac. */
223 real shake_tol; /* tolerance for shake */
224 int efep; /* free energy interpolation no/yes */
225 double init_lambda; /* initial value for perturbation variable */
226 double delta_lambda; /* change of lambda per time step (1/dt) */
227 int n_flambda; /* The number of foreign lambda points */
228 double *flambda; /* The foreign lambda values */
229 real sc_alpha; /* free energy soft-core parameter */
230 int sc_power; /* lambda power for soft-core interactions */
231 real sc_sigma; /* free energy soft-core sigma when c6 or c12=0 */
232 int nstdhdl; /* The frequency for writing to dhdl.xvg */
233 int dh_table_size; /* The maximum table size for the dH table */
234 double dh_table_spacing; /* The spacing for the dH table */
235 int eDisre; /* Type of distance restraining */
236 real dr_fc; /* force constant for ta_disre */
237 int eDisreWeighting; /* type of weighting of pairs in one restraints */
238 bool bDisreMixed; /* Use comb of time averaged and instan. viol's */
239 int nstdisreout; /* frequency of writing pair distances to enx */
240 real dr_tau; /* time constant for memory function in disres */
241 real orires_fc; /* force constant for orientational restraints */
242 real orires_tau; /* time constant for memory function in orires */
243 int nstorireout; /* frequency of writing tr(SD) to enx */
244 real dihre_fc; /* force constant for dihedral restraints */
245 real em_stepsize; /* The stepsize for updating */
246 real em_tol; /* The tolerance */
247 int niter; /* Number of iterations for convergence of */
248 /* steepest descent in relax_shells */
249 real fc_stepsize; /* Stepsize for directional minimization */
250 /* in relax_shells */
251 int nstcgsteep; /* number of steps after which a steepest */
252 /* descents step is done while doing cg */
253 int nbfgscorr; /* Number of corrections to the hessian to keep */
254 int eConstrAlg; /* Type of constraint algorithm */
255 int nProjOrder; /* Order of the LINCS Projection Algorithm */
256 real LincsWarnAngle; /* If bond rotates more than %g degrees, warn */
257 int nLincsIter; /* Number of iterations in the final Lincs step */
258 bool bShakeSOR; /* Use successive overrelaxation for shake */
259 real bd_fric; /* Friction coefficient for BD (amu/ps) */
260 int ld_seed; /* Random seed for SD and BD */
261 int nwall; /* The number of walls */
262 int wall_type; /* The type of walls */
263 real wall_r_linpot; /* The potentail is linear for r<=wall_r_linpot */
264 int wall_atomtype[2];/* The atom type for walls */
265 real wall_density[2]; /* Number density for walls */
266 real wall_ewald_zfac; /* Scaling factor for the box for Ewald */
267 int ePull; /* Type of pulling: no, umbrella or constraint */
268 t_pull *pull; /* The data for center of mass pulling */
269 real cos_accel; /* Acceleration for viscosity calculation */
270 tensor deform; /* Triclinic deformation velocities (nm/ps) */
271 int userint1; /* User determined parameters */
272 int userint2;
273 int userint3;
274 int userint4;
275 real userreal1;
276 real userreal2;
277 real userreal3;
278 real userreal4;
279 t_grpopts opts; /* Group options */
280 t_cosines ex[DIM]; /* Electric field stuff (spatial part) */
281 t_cosines et[DIM]; /* Electric field stuff (time part) */
282 bool bQMMM; /* QM/MM calculation */
283 int QMconstraints; /* constraints on QM bonds */
284 int QMMMscheme; /* Scheme: ONIOM or normal */
285 real scalefactor; /* factor for scaling the MM charges in QM calc.*/
286 } t_inputrec;
288 #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)
290 #define DYNAMIC_BOX(ir) ((ir).epc!=epcNO || (ir).eI==eiTPI || DEFORM(ir))
292 #define PRESERVE_SHAPE(ir) ((ir).epc != epcNO && (ir).deform[XX][XX] == 0 && ((ir).epct == epctISOTROPIC || (ir).epct == epctSEMIISOTROPIC))
294 #define NEED_MUTOT(ir) (((ir).coulombtype==eelEWALD || EEL_PME((ir).coulombtype)) && ((ir).ewald_geometry==eewg3DC || (ir).epsilon_surface!=0))
296 #define IR_TWINRANGE(ir) ((ir).rlist > 0 && ((ir).rlistlong == 0 || (ir).rlistlong > (ir).rlist))
298 #define IR_ELEC_FIELD(ir) ((ir).ex[XX].n > 0 || (ir).ex[YY].n > 0 || (ir).ex[ZZ].n > 0)
300 #define IR_EXCL_FORCES(ir) (EEL_FULL((ir).coulombtype) || (EEL_RF((ir).coulombtype) && (ir).coulombtype != eelRF_NEC) || (ir).implicit_solvent != eisNO)
301 /* use pointer definitions of ir here, since that's what's usually used in the code */
302 #define IR_NVT_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && ((ir)->etc == etcNOSEHOOVER))
304 #define IR_NPT_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && ((ir)->epc == epcMTTK))
306 #ifdef __cplusplus
308 #endif
311 #endif