Tagged files with gromacs 3.0 header and added some license info
[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 3.0
11  * 
12  * Copyright (c) 1991-2001
13  * BIOSON Research Institute, Dept. of Biophysical Chemistry
14  * University of Groningen, The Netherlands
15  * 
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
32  * 
33  * And Hey:
34  * Good ROcking Metal Altar for Chronical Sinners
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 typedef struct {
41   int  n;               /* Number of terms                              */
42   real *a;              /* Coeffients (V / nm )                         */
43   real *phi;            /* Phase angles                                 */
44 } t_cosines;
45
46 typedef struct {
47   int     ngtc;         /* # T-Coupl groups                             */
48   int     ngacc;        /* # Accelerate groups                          */
49   int     ngfrz;        /* # Freeze groups                              */
50   int     ngener;       /* # Ener groups                                */
51   real    *nrdf;        /* Nr of degrees of freedom in a group          */
52   real    *ref_t;       /* Coupling temperature per group               */
53   real    *tau_t;       /* Tau coupling time                            */
54   rvec    *acc;         /* Acceleration per group                       */
55   ivec    *nFreeze;     /* Freeze the group in each direction ?         */
56   bool    *eg_excl;     /* Exclusions of energy group pairs             */
57 } t_grpopts;
58
59 typedef struct {
60   int  eI;              /* Integration method                           */
61   int  nsteps;          /* number of steps to be taken                  */
62   int  ns_type;         /* which ns method should we use?               */
63   int  nstlist;         /* number of steps before pairlist is generated */
64   int  ndelta;          /* number of cells per rlong                    */
65   bool bDomDecomp;      /* Should we do domain decomposition?           */
66   int  decomp_dir;      /* Direction of decomposition (may not be opt.) */
67   int  nstcomm;         /* number of steps after which center of mass   */
68                         /* motion is removed                            */
69   int nstlog;           /* number of steps after which print to logfile */
70   int nstxout;          /* number of steps after which X is output      */
71   int nstvout;          /* id. for V                                    */
72   int nstfout;          /* id. for F                                    */
73   int nstenergy;        /* number of steps after which energies printed */
74   int nstxtcout;        /* id. for compressed trj (.xtc)                */
75   real init_t;          /* initial time (ps)                            */
76   real delta_t;         /* time step (ps)                               */
77   real xtcprec;         /* precision of xtc file                        */
78   int  nkx,nky,nkz;     /* number of k vectors in each spatial dimension*/
79                         /* for fourier methods for long range electrost.*/
80   int  pme_order;       /* interpolation order for PME                  */
81   real ewald_rtol;      /* Real space tolerance for Ewald, determines   */
82                         /* the real/reciprocal space relative weight    */
83   bool epsilon_surface; /* Epsilon for PME dipole correction            */
84   bool bOptFFT;         /* optimize the fft plan at start               */
85   int  ePBC;            /* Type of periodic boundary conditions         */
86   bool bUncStart;       /* Do not constrain the start configuration     */
87   int  etc;             /* temperature coupling                         */
88   int  epc;             /* pressure coupling                            */
89   int  epct;            /* pressure coupling type                       */
90   real tau_p;           /* pressure coupling time (ps)                  */
91   tensor ref_p;         /* reference pressure (kJ/(mol nm^3))           */
92   tensor compress;      /* compressability ((mol nm^3)/kJ)              */
93   bool bSimAnn;         /* simulated annealing (SA)                     */
94   real zero_temp_time;  /* time at which temp becomes zero in sim. ann. */
95   real rlist;           /* short range pairlist cut-off (nm)            */
96   int  coulombtype;     /* Type of electrostatics treatment             */
97   real rcoulomb_switch; /* Coulomb switch range start (nm)              */
98   real rcoulomb;        /* Coulomb cutoff (nm)                          */
99   int  vdwtype;         /* Type of Van der Waals treatment              */
100   real rvdw_switch;     /* Van der Waals switch range start (nm)        */
101   real rvdw;            /* Van der Waals cutoff (nm)                    */
102   real epsilon_r;       /* relative dielectric constant                 */
103   int  eDispCorr;       /* Perform Long range dispersion corrections    */
104   real shake_tol;       /* tolerance for shake                          */
105   real fudgeQQ;         /* Id. for 1-4 coulomb interactions             */
106   int  efep;            /* free energy interpolation no/yes             */
107   real init_lambda;     /* initial value for perturbation variable      */
108   real delta_lambda;    /* change of lambda per time step (1/dt)        */
109   real sc_alpha;        /* free energy soft-core parameter              */
110   real sc_sigma;        /* free energy soft-core sigma when c6 or c12=0 */
111   real dr_fc;           /* force constant for ta_disre                  */
112   int  eDisreWeighting; /* type of weighting of pairs in one restraints */
113   bool bDisreMixed;     /* Use comb of time averaged and instan. viol's */
114   int  nstdisreout;     /* frequency of writing pair distances to enx   */ 
115   real dr_tau;          /* time constant for memory function in disres  */
116   real em_stepsize;     /* The stepsize for updating                    */
117   real em_tol;          /* The tolerance                                */
118   int  niter;           /* Number of iterations for convergence of      */
119                         /* steepest descent in relax_shells             */
120   int  nstcgsteep;      /* number of steps after which a steepest       */
121                         /* descents step is done while doing cg         */
122   int  eConstrAlg;      /* Type of constraint algorithm                 */
123   int  nProjOrder;      /* Order of the LINCS Projection Algorithm      */
124   real LincsWarnAngle;  /* If bond rotates more than %g degrees, warn   */
125   real bd_temp;         /* Temperature for Brownian Dynamics (BD)       */
126   real bd_fric;         /* Friction coefficient for BD (amu / ps)       */
127   int  ld_seed;         /* Random seed for SD and BD                    */
128   real cos_accel;       /* Acceleration for viscosity calculation       */
129   int  userint1;        /* User determined parameters                   */
130   int  userint2;
131   int  userint3;
132   int  userint4;
133   real userreal1;
134   real userreal2;
135   real userreal3;
136   real userreal4;
137   t_grpopts opts;       /* Group options                                */
138   t_cosines ex[DIM];    /* Electric field stuff (spatial part)          */
139   t_cosines et[DIM];    /* Electric field stuff (time part)             */
140 } t_inputrec;
141