992cf3a47bc29c8b34a61c4a3055491f6c9d17a4
[alexxy/gromacs.git] / include / types / inputrec.h
1 /*
2  * $Id$
3  * 
4  *       This source code is part of
5  * 
6  *        G   R   O   M   A   C   S
7  * 
8  * GROningen MAchine for Chemical Simulations
9  * 
10  *               VERSION 2.0
11  * 
12  * Copyright (c) 1991-1999
13  * BIOSON Research Institute, Dept. of Biophysical Chemistry
14  * University of Groningen, The Netherlands
15  * 
16  * Please refer to:
17  * GROMACS: A message-passing parallel molecular dynamics implementation
18  * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19  * Comp. Phys. Comm. 91, 43-56 (1995)
20  * 
21  * Also check out our WWW page:
22  * http://md.chem.rug.nl/~gmx
23  * or e-mail to:
24  * gromacs@chem.rug.nl
25  * 
26  * And Hey:
27  * Green Red Orange Magenta Azure Cyan Skyblue
28  */
29
30 #ifdef HAVE_CONFIG_H
31 #include <config.h>
32 #endif
33
34 typedef struct {
35   int  n;               /* Number of terms                              */
36   real *a;              /* Coeffients (V / nm )                         */
37   real *phi;            /* Phase angles                                 */
38 } t_cosines;
39
40 typedef struct {
41   int     ngtc;         /* # T-Coupl groups                             */
42   int     ngacc;        /* # Accelerate groups                          */
43   int     ngfrz;        /* # Freeze groups                              */
44   int     ngener;       /* # Ener groups                                */
45   real    *nrdf;        /* Nr of degrees of freedom in a group          */
46   real    *ref_t;       /* Coupling temperature per group               */
47   real    *tau_t;       /* Tau coupling time                            */
48   rvec    *acc;         /* Acceleration per group                       */
49   ivec    *nFreeze;     /* Freeze the group in each direction ?         */
50   bool    *eg_excl;     /* Exclusions of energy group pairs             */
51 } t_grpopts;
52
53 typedef struct {
54   int  eI;              /* Integration method                           */
55   int  nsteps;          /* number of steps to be taken                  */
56   int  ns_type;         /* which ns method should we use?               */
57   int  nstlist;         /* number of steps before pairlist is generated */
58   int  ndelta;          /* number of cells per rlong                    */
59   bool bDomDecomp;      /* Should we do domain decomposition?           */
60   int  decomp_dir;      /* Direction of decomposition (may not be opt.) */
61   int  nstcomm;         /* number of steps after which center of mass   */
62                         /* motion is removed                            */
63   int nstlog;           /* number of steps after which print to logfile */
64   int nstxout;          /* number of steps after which X is output      */
65   int nstvout;          /* id. for V                                    */
66   int nstfout;          /* id. for F                                    */
67   int nstenergy;        /* number of steps after which energies printed */
68   int nstxtcout;        /* id. for compressed trj (.xtc)                */
69   real init_t;          /* initial time (ps)                            */
70   real delta_t;         /* time step (ps)                               */
71   real xtcprec;         /* precision of xtc file                        */
72   int  nkx,nky,nkz;     /* number of k vectors in each spatial dimension*/
73                         /* for fourier methods for long range electrost.*/
74   int  pme_order;       /* interpolation order for PME                  */
75   real ewald_rtol;      /* Real space tolerance for Ewald, determines   */
76                         /* the real/reciprocal space relative weight    */
77   bool epsilon_surface; /* Epsilon for PME dipole correction            */
78   bool bOptFFT;         /* optimize the fft plan at start               */
79   int  ePBC;            /* Type of periodic boundary conditions         */
80   bool bUncStart;       /* Do not constrain the start configuration     */
81   int  etc;             /* temperature coupling                         */
82   int  epc;             /* pressure coupling                            */
83   int  epct;            /* pressure coupling type                       */
84   real tau_p;           /* pressure coupling time (ps)                  */
85   tensor ref_p;         /* reference pressure (kJ/(mol nm^3))           */
86   tensor compress;      /* compressability ((mol nm^3)/kJ)              */
87   bool bSimAnn;         /* simulated annealing (SA)                     */
88   real zero_temp_time;  /* time at which temp becomes zero in sim. ann. */
89   real rlist;           /* short range pairlist cut-off (nm)            */
90   int  coulombtype;     /* Type of electrostatics treatment             */
91   real rcoulomb_switch; /* Coulomb switch range start (nm)              */
92   real rcoulomb;        /* Coulomb cutoff (nm)                          */
93   int  vdwtype;         /* Type of Van der Waals treatment              */
94   real rvdw_switch;     /* Van der Waals switch range start (nm)        */
95   real rvdw;            /* Van der Waals cutoff (nm)                    */
96   real epsilon_r;       /* relative dielectric constant                 */
97   int  eDispCorr;       /* Perform Long range dispersion corrections    */
98   real shake_tol;       /* tolerance for shake                          */
99   real fudgeQQ;         /* Id. for 1-4 coulomb interactions             */
100   int  efep;            /* free energy interpolation no/yes             */
101   real init_lambda;     /* initial value for perturbation variable      */
102   real delta_lambda;    /* change of lambda per time step (1/dt)        */
103   real sc_alpha;        /* free energy soft-core parameter              */
104   real sc_sigma;        /* free energy soft-core sigma when c6 or c12=0 */
105   real dr_fc;           /* force constant for ta_disre                  */
106   int  eDisreWeighting; /* type of weighting of pairs in one restraints */
107   bool bDisreMixed;     /* Use comb of time averaged and instan. viol's */
108   int  nstdisreout;     /* frequency of writing pair distances to enx   */ 
109   real dr_tau;          /* time constant for memory function in disres  */
110   real em_stepsize;     /* The stepsize for updating                    */
111   real em_tol;          /* The tolerance                                */
112   int  niter;           /* Number of iterations for convergence of      */
113                         /* steepest descent in relax_shells             */
114   int  nstcgsteep;      /* number of steps after which a steepest       */
115                         /* descents step is done while doing cg         */
116   int  eConstrAlg;      /* Type of constraint algorithm                 */
117   int  nProjOrder;      /* Order of the LINCS Projection Algorithm      */
118   real LincsWarnAngle;  /* If bond rotates more than %g degrees, warn   */
119   real bd_temp;         /* Temperature for Brownian Dynamics (BD)       */
120   real bd_fric;         /* Friction coefficient for BD (amu / ps)       */
121   int  ld_seed;         /* Random seed for SD and BD                    */
122   real cos_accel;       /* Acceleration for viscosity calculation       */
123   int  userint1;        /* User determined parameters                   */
124   int  userint2;
125   int  userint3;
126   int  userint4;
127   real userreal1;
128   real userreal2;
129   real userreal3;
130   real userreal4;
131   t_grpopts opts;       /* Group options                                */
132   t_cosines ex[DIM];    /* Electric field stuff (spatial part)          */
133   t_cosines et[DIM];    /* Electric field stuff (time part)             */
134 } t_inputrec;
135