changed reading hint
[gromacs/adressmacs.git] / include / types / forcerec.h
blob9ac5a1c3ad879e961a5ba7f18ff72f16188c7339
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
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
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)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * Green Red Orange Magenta Azure Cyan Skyblue
30 enum { eNL_VDW, eNL_QQ, eNL_FREE, eNL_VDW_WAT, eNL_QQ_WAT, eNL_NR };
32 typedef struct {
33 /* Cut-Off stuff */
34 int eBox;
35 real rlist,rlistlong;
37 /* Dielectric constant resp. multiplication factor for charges */
38 real zsquare,temp;
39 real epsilon_r,epsfac;
41 /* Constants for reaction fields */
42 bool bRF;
43 real kappa,k_rf,c_rf;
45 /* Constant for long range dispersion correction (average dispersion) */
46 real avcsix;
48 /* Fudge factors */
49 real fudgeQQ;
51 /* Table stuff */
52 bool bTab;
53 real rtab;
54 int ntab;
55 real tabscale;
56 real *VFtab;
57 real *VFtab14; /* special 1-4 tables to use when the ordinary
58 * tabulated coulomb interactions are LR shifted.
59 * not a pretty solution to double LJ, but it works.
60 * Erik 990903
63 /* PPPM & Shifting stuff */
64 real rcoulomb_switch,rcoulomb;
65 real *phi;
67 /* VdW stuff */
68 real rvdw_switch,rvdw;
69 real bham_b_max;
70 real tabscale_exp;
72 /* Free energy ? */
73 bool bPert;
75 /* NS Stuff */
76 int eeltype;
77 int vdwtype;
78 int cg0,hcg;
79 int ndelta;
80 bool bWaterOpt;
81 int nWater,nWatMol;
82 int Dimension;
83 bool bGrid,bDomDecomp;
84 rvec *cg_cm;
85 rvec *shift_vec;
87 /* The actual neighbor lists, short and long range, see enum above
88 * for definition of neighborlist indices.
90 t_nblist nlist_sr[eNL_NR];
91 t_nblist nlist_lr[eNL_NR];
93 /* This mask array of length nn determines whether or not this bit of the
94 * neighbourlists should be computed. Usually all these are true of course,
95 * but not when shells are used. During minimisation all the forces that
96 * include shells are done, then after minimsation is converged the remaining
97 * forces are computed.
99 /* bool *bMask; */
101 /* Twin Range stuff. */
102 bool bTwinRange;
103 int nlr;
104 rvec *flr;
105 rvec *fshift_lr;
107 /* PME/Ewald stuff */
108 bool bEwald;
109 real ewaldcoeff;
111 /* Virial Stuff */
112 rvec *fshift;
114 /* Free energy stuff */
115 int nmol;
116 atom_id *mol_nr;
117 real *mol_epot;
118 int nstcalc;
120 /* Non bonded Parameter lists */
121 int ntype; /* Number of atom types */
122 bool bBHAM;
123 real *nbfp;
124 } t_forcerec;
126 #define C6(nbfp,ntp,ai,aj) (nbfp)[2*((ntp)*(ai)+(aj))]
127 #define C12(nbfp,ntp,ai,aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
128 #define BHAMC(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))]
129 #define BHAMA(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
130 #define BHAMB(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]