3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * GRoups of Organic Molecules in ACtion for Science
40 #include "../sysstuff.h"
48 int n; /* Number of terms */
49 real *a; /* Coeffients (V / nm ) */
50 real *phi; /* Phase angles */
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) */
60 #define EGP_EXCL (1<<0)
61 #define EGP_TABLE (1<<1)
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 */
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 */
92 int *SAsteps; /* in how many steps SA goes from 1-1 to 0.5-0.5*/
97 enum { epgrppbcNONE, epgrppbcREFAT, epgrppbcCOS };
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 */
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 */
132 int nstdhdl; /* The frequency for calculating dhdl */
133 double init_lambda; /* fractional value of lambda (usually will use
134 init_fep_state, this will only be for slow growth,
135 and for legacy free energy code. Only has a
136 valid value if positive) */
137 int init_fep_state; /* the initial number of the state */
138 double delta_lambda; /* change of lambda per time step (fraction of (0.1) */
139 gmx_bool bPrintEnergy; /* Whether to print the energy in the dhdl */
140 int n_lambda; /* The number of foreign lambda points */
141 double **all_lambda; /* The array of all lambda values */
142 int lambda_neighbors; /* The number of neighboring lambda states to
143 calculate the energy for in up and down directions
145 int lambda_start_n; /* The first lambda to calculate energies for */
146 int lambda_stop_n; /* The last lambda +1 to calculate energies for */
147 real sc_alpha; /* free energy soft-core parameter */
148 int sc_power; /* lambda power for soft-core interactions */
149 real sc_r_power; /* r power for soft-core interactions */
150 real sc_sigma; /* free energy soft-core sigma when c6 or c12=0 */
151 real sc_sigma_min; /* free energy soft-core sigma for ????? */
152 gmx_bool bScCoul; /* use softcore for the coulomb portion as well (default FALSE) */
153 gmx_bool separate_dvdl[efptNR]; /* whether to print the dvdl term associated with
154 this term; if it is not specified as separate,
155 it is lumped with the FEP term */
156 int separate_dhdl_file; /* whether to write a separate dhdl.xvg file
157 note: NOT a gmx_bool, but an enum */
158 int dhdl_derivatives; /* whether to calculate+write dhdl derivatives
159 note: NOT a gmx_bool, but an enum */
160 int dh_hist_size; /* The maximum table size for the dH histogram */
161 double dh_hist_spacing; /* The spacing for the dH histogram */
165 int nstexpanded; /* The frequency of expanded ensemble state changes */
166 int elamstats; /* which type of move updating do we use for lambda monte carlo (or no for none) */
167 int elmcmove; /* what move set will be we using for state space moves */
168 int elmceq; /* the method we use to decide of we have equilibrated the weights */
169 int equil_n_at_lam; /* the minumum number of samples at each lambda for deciding whether we have reached a minimum */
170 real equil_wl_delta; /* WL delta at which we stop equilibrating weights */
171 real equil_ratio; /* use the ratio of weights (ratio of minimum to maximum) to decide when to stop equilibrating */
172 int equil_steps; /* after equil_steps steps we stop equilibrating the weights */
173 int equil_samples; /* after equil_samples total samples (steps/nstfep), we stop equilibrating the weights */
174 int lmc_seed; /* random number seed for lambda mc switches */
175 gmx_bool minvar; /* whether to use minumum variance weighting */
176 int minvarmin; /* the number of samples needed before kicking into minvar routine */
177 real minvar_const; /* the offset for the variance in MinVar */
178 int c_range; /* range of cvalues used for BAR */
179 gmx_bool bSymmetrizedTMatrix; /* whether to print symmetrized matrices */
180 int nstTij; /* How frequently to print the transition matrices */
181 int lmc_repeats; /* number of repetitions in the MC lambda jumps */ /*MRS -- VERIFY THIS */
182 int lmc_forced_nstart; /* minimum number of samples for each state before free sampling */ /* MRS -- VERIFY THIS! */
183 int gibbsdeltalam; /* distance in lambda space for the gibbs interval */
184 real wl_scale; /* scaling factor for wang-landau */
185 real wl_ratio; /* ratio between largest and smallest number for freezing the weights */
186 real init_wl_delta; /* starting delta for wang-landau */
187 gmx_bool bWLoneovert; /* use one over t convergece for wang-landau when the delta get sufficiently small */
188 gmx_bool bInit_weights; /* did we initialize the weights? */
189 real mc_temp; /* To override the main temperature, or define it if it's not defined */
190 real *init_lambda_weights; /* user-specified initial weights to start with */
194 int ngrp; /* number of groups */
195 int eGeom; /* pull geometry */
196 ivec dim; /* used to select components for constraint */
197 real cyl_r1; /* radius of cylinder for dynamic COM */
198 real cyl_r0; /* radius of cylinder including switch length */
199 real constr_tol; /* absolute tolerance for constraints in (nm) */
200 int nstxout; /* Output frequency for pull x */
201 int nstfout; /* Output frequency for pull f */
202 int ePBC; /* the boundary conditions */
203 int npbcdim; /* do pbc in dims 0 <= dim < npbcdim */
204 gmx_bool bRefAt; /* do we need reference atoms for a group COM ? */
205 int cosdim; /* dimension for cosine weighting, -1 if none */
206 gmx_bool bVirial; /* do we need to add the pull virial? */
207 t_pullgrp *grp; /* groups to pull/restrain/etc/ */
208 t_pullgrp *dyna; /* dynamic groups for use with local constraints */
209 rvec *rbuf; /* COM calculation buffer */
210 dvec *dbuf; /* COM calculation buffer */
211 double *dbuf_cyl; /* cylinder ref. groups COM calculation buffer */
213 FILE *out_x; /* output file for pull data */
214 FILE *out_f; /* output file for pull data */
218 /* Abstract types for enforced rotation only defined in pull_rotation.c */
219 typedef struct gmx_enfrot *gmx_enfrot_t;
220 typedef struct gmx_enfrotgrp *gmx_enfrotgrp_t;
223 int eType; /* Rotation type for this group */
224 int bMassW; /* Use mass-weighed positions? */
225 int nat; /* Number of atoms in the group */
226 atom_id *ind; /* The global atoms numbers */
227 rvec *x_ref; /* The reference positions */
228 rvec vec; /* The normalized rotation vector */
229 real rate; /* Rate of rotation (degree/ps) */
230 real k; /* Force constant (kJ/(mol nm^2) */
231 rvec pivot; /* Pivot point of rotation axis (nm) */
232 int eFittype; /* Type of fit to determine actual group angle */
233 int PotAngle_nstep; /* Number of angles around the reference angle
234 for which the rotation potential is also
235 evaluated (for fit type 'potential' only) */
236 real PotAngle_step; /* Distance between two angles in degrees (for
237 fit type 'potential' only) */
238 real slab_dist; /* Slab distance (nm) */
239 real min_gaussian; /* Minimum value the gaussian must have so that
240 the force is actually evaluated */
241 real eps; /* Additive constant for radial motion2 and
242 flexible2 potentials (nm^2) */
243 gmx_enfrotgrp_t enfrotgrp; /* Stores non-inputrec rotation data per group */
247 int ngrp; /* Number of rotation groups */
248 int nstrout; /* Output frequency for main rotation outfile */
249 int nstsout; /* Output frequency for per-slab data */
250 t_rotgrp *grp; /* Groups to rotate */
251 gmx_enfrot_t enfrot; /* Stores non-inputrec enforced rotation data */
256 int type; /* type of AdResS simulation */
257 gmx_bool bnew_wf; /* enable new AdResS weighting function */
258 gmx_bool bchempot_dx; /*true:interaction table format input is F=-dmu/dx false: dmu_dwp */
259 gmx_bool btf_full_box; /* true: appy therm force everywhere in the box according to table false: only in hybrid region */
260 real const_wf; /* value of weighting function for eAdressConst */
261 real ex_width; /* center of the explicit zone */
262 real hy_width; /* width of the hybrid zone */
263 int icor; /* type of interface correction */
264 int site; /* AdResS CG site location */
265 rvec refs; /* Coordinates for AdResS reference */
266 real ex_forcecap; /* in the hybrid zone, cap forces large then this to adress_ex_forcecap */
267 gmx_bool do_hybridpairs; /* If true pair interaction forces are also scaled in an adress way*/
269 int * tf_table_index; /* contains mapping of energy group index -> i-th adress tf table*/
276 int eI; /* Integration method */
277 gmx_large_int_t nsteps; /* number of steps to be taken */
278 int simulation_part; /* Used in checkpointing to separate chunks */
279 gmx_large_int_t init_step; /* start at a stepcount >0 (used w. tpbconv) */
280 int nstcalcenergy; /* frequency of energy calc. and T/P coupl. upd. */
281 int cutoff_scheme; /* group or verlet cutoffs */
282 int ns_type; /* which ns method should we use? */
283 int nstlist; /* number of steps before pairlist is generated */
284 int ndelta; /* number of cells per rlong */
285 int nstcomm; /* number of steps after which center of mass */
286 /* motion is removed */
287 int comm_mode; /* Center of mass motion removal algorithm */
288 int nstcheckpoint; /* checkpointing frequency */
289 int nstlog; /* number of steps after which print to logfile */
290 int nstxout; /* number of steps after which X is output */
291 int nstvout; /* id. for V */
292 int nstfout; /* id. for F */
293 int nstenergy; /* number of steps after which energies printed */
294 int nstxtcout; /* id. for compressed trj (.xtc) */
295 double init_t; /* initial time (ps) */
296 double delta_t; /* time step (ps) */
297 real xtcprec; /* precision of xtc file */
298 real fourier_spacing; /* requested fourier_spacing, when nk? not set */
299 int nkx,nky,nkz; /* number of k vectors in each spatial dimension*/
300 /* for fourier methods for long range electrost.*/
301 int pme_order; /* interpolation order for PME */
302 real ewald_rtol; /* Real space tolerance for Ewald, determines */
303 /* the real/reciprocal space relative weight */
304 int ewald_geometry; /* normal/3d ewald, or pseudo-2d LR corrections */
305 real epsilon_surface; /* Epsilon for PME dipole correction */
306 gmx_bool bOptFFT; /* optimize the fft plan at start */
307 int ePBC; /* Type of periodic boundary conditions */
308 int bPeriodicMols; /* Periodic molecules */
309 gmx_bool bContinuation; /* Continuation run: starting state is correct */
310 int etc; /* temperature coupling */
311 int nsttcouple; /* interval in steps for temperature coupling */
312 gmx_bool bPrintNHChains; /* whether to print nose-hoover chains */
313 int epc; /* pressure coupling */
314 int epct; /* pressure coupling type */
315 int nstpcouple; /* interval in steps for pressure coupling */
316 real tau_p; /* pressure coupling time (ps) */
317 tensor ref_p; /* reference pressure (kJ/(mol nm^3)) */
318 tensor compress; /* compressability ((mol nm^3)/kJ) */
319 int refcoord_scaling;/* How to scale absolute reference coordinates */
320 rvec posres_com; /* The COM of the posres atoms */
321 rvec posres_comB; /* The B-state COM of the posres atoms */
322 int andersen_seed; /* Random seed for Andersen thermostat (obsolete) */
323 real verletbuf_drift; /* Max. drift (kJ/mol/ps/atom) for list buffer */
324 real rlist; /* short range pairlist cut-off (nm) */
325 real rlistlong; /* long range pairlist cut-off (nm) */
326 int nstcalclr; /* Frequency of evaluating direct space long-range interactions */
327 real rtpi; /* Radius for test particle insertion */
328 int coulombtype; /* Type of electrostatics treatment */
329 int coulomb_modifier; /* Modify the Coulomb interaction */
330 real rcoulomb_switch; /* Coulomb switch range start (nm) */
331 real rcoulomb; /* Coulomb cutoff (nm) */
332 real epsilon_r; /* relative dielectric constant */
333 real epsilon_rf; /* relative dielectric constant of the RF */
334 int implicit_solvent;/* No (=explicit water), or GBSA solvent models */
335 int gb_algorithm; /* Algorithm to use for calculation Born radii */
336 int nstgbradii; /* Frequency of updating Generalized Born radii */
337 real rgbradii; /* Cutoff for GB radii calculation */
338 real gb_saltconc; /* Salt concentration (M) for GBSA models */
339 real gb_epsilon_solvent; /* dielectric coeff. of implicit solvent */
340 real gb_obc_alpha; /* 1st scaling factor for Bashford-Case GB */
341 real gb_obc_beta; /* 2nd scaling factor for Bashford-Case GB */
342 real gb_obc_gamma; /* 3rd scaling factor for Bashford-Case GB */
343 real gb_dielectric_offset; /* Dielectric offset for Still/HCT/OBC */
344 int sa_algorithm; /* Algorithm for SA part of GBSA */
345 real sa_surface_tension; /* Energy factor for SA part of GBSA */
346 int vdwtype; /* Type of Van der Waals treatment */
347 int vdw_modifier; /* Modify the VdW interaction */
348 real rvdw_switch; /* Van der Waals switch range start (nm) */
349 real rvdw; /* Van der Waals cutoff (nm) */
350 int eDispCorr; /* Perform Long range dispersion corrections */
351 real tabext; /* Extension of the table beyond the cut-off, *
352 * as well as the table length for 1-4 interac. */
353 real shake_tol; /* tolerance for shake */
354 int efep; /* free energy calculations */
355 t_lambda *fepvals; /* Data for the FEP state */
356 gmx_bool bSimTemp; /* Whether to do simulated tempering */
357 t_simtemp *simtempvals;/* Variables for simulated tempering */
358 gmx_bool bExpanded; /* Whether expanded ensembles are used */
359 t_expanded *expandedvals; /* Expanded ensemble parameters */
360 int eDisre; /* Type of distance restraining */
361 real dr_fc; /* force constant for ta_disre */
362 int eDisreWeighting; /* type of weighting of pairs in one restraints */
363 gmx_bool bDisreMixed; /* Use comb of time averaged and instan. viol's */
364 int nstdisreout; /* frequency of writing pair distances to enx */
365 real dr_tau; /* time constant for memory function in disres */
366 real orires_fc; /* force constant for orientational restraints */
367 real orires_tau; /* time constant for memory function in orires */
368 int nstorireout; /* frequency of writing tr(SD) to enx */
369 real dihre_fc; /* force constant for dihedral restraints (obsolete) */
370 real em_stepsize; /* The stepsize for updating */
371 real em_tol; /* The tolerance */
372 int niter; /* Number of iterations for convergence of */
373 /* steepest descent in relax_shells */
374 real fc_stepsize; /* Stepsize for directional minimization */
375 /* in relax_shells */
376 int nstcgsteep; /* number of steps after which a steepest */
377 /* descents step is done while doing cg */
378 int nbfgscorr; /* Number of corrections to the hessian to keep */
379 int eConstrAlg; /* Type of constraint algorithm */
380 int nProjOrder; /* Order of the LINCS Projection Algorithm */
381 real LincsWarnAngle; /* If bond rotates more than %g degrees, warn */
382 int nLincsIter; /* Number of iterations in the final Lincs step */
383 gmx_bool bShakeSOR; /* Use successive overrelaxation for shake */
384 real bd_fric; /* Friction coefficient for BD (amu/ps) */
385 int ld_seed; /* Random seed for SD and BD */
386 int nwall; /* The number of walls */
387 int wall_type; /* The type of walls */
388 real wall_r_linpot; /* The potentail is linear for r<=wall_r_linpot */
389 int wall_atomtype[2];/* The atom type for walls */
390 real wall_density[2]; /* Number density for walls */
391 real wall_ewald_zfac; /* Scaling factor for the box for Ewald */
392 int ePull; /* Type of pulling: no, umbrella or constraint */
393 t_pull *pull; /* The data for center of mass pulling */
394 gmx_bool bRot; /* Calculate enforced rotation potential(s)? */
395 t_rot *rot; /* The data for enforced rotation potentials */
396 real cos_accel; /* Acceleration for viscosity calculation */
397 tensor deform; /* Triclinic deformation velocities (nm/ps) */
398 int userint1; /* User determined parameters */
406 t_grpopts opts; /* Group options */
407 t_cosines ex[DIM]; /* Electric field stuff (spatial part) */
408 t_cosines et[DIM]; /* Electric field stuff (time part) */
409 gmx_bool bQMMM; /* QM/MM calculation */
410 int QMconstraints; /* constraints on QM bonds */
411 int QMMMscheme; /* Scheme: ONIOM or normal */
412 real scalefactor; /* factor for scaling the MM charges in QM calc.*/
413 /* parameter needed for AdResS simulation */
414 gmx_bool bAdress; /* Is AdResS enabled ? */
415 t_adress *adress; /* The data for adress simulations */
418 #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)
420 #define DYNAMIC_BOX(ir) ((ir).epc!=epcNO || (ir).eI==eiTPI || DEFORM(ir))
422 #define PRESERVE_SHAPE(ir) ((ir).epc != epcNO && (ir).deform[XX][XX] == 0 && ((ir).epct == epctISOTROPIC || (ir).epct == epctSEMIISOTROPIC))
424 #define NEED_MUTOT(ir) (((ir).coulombtype==eelEWALD || EEL_PME((ir).coulombtype)) && ((ir).ewald_geometry==eewg3DC || (ir).epsilon_surface!=0))
426 #define IR_TWINRANGE(ir) ((ir).rlist > 0 && ((ir).rlistlong == 0 || (ir).rlistlong > (ir).rlist))
428 #define IR_ELEC_FIELD(ir) ((ir).ex[XX].n > 0 || (ir).ex[YY].n > 0 || (ir).ex[ZZ].n > 0)
430 #define IR_EXCL_FORCES(ir) (EEL_FULL((ir).coulombtype) || (EEL_RF((ir).coulombtype) && (ir).coulombtype != eelRF_NEC) || (ir).implicit_solvent != eisNO)
431 /* use pointer definitions of ir here, since that's what's usually used in the code */
432 #define IR_NPT_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && (((ir)->epc == epcMTTK) && ((ir)->etc == etcNOSEHOOVER)))
434 #define IR_NVT_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && ((!((ir)->epc == epcMTTK)) && ((ir)->etc == etcNOSEHOOVER)))
436 #define IR_NPH_TROTTER(ir) ((((ir)->eI == eiVV) || ((ir)->eI == eiVVAK)) && (((ir)->epc == epcMTTK) && (!(((ir)->etc == etcNOSEHOOVER)))))