added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / include / types / idef.h
blob7d44e4ccb2a881627feb6a51ec03d20fb22ba382
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
37 #ifndef _idef_h
38 #define _idef_h
40 #include "simple.h"
42 #ifdef __cplusplus
43 extern "C" {
44 #endif
47 /* check kernel/toppush.c when you change these numbers */
48 #define MAXATOMLIST 6
49 #define MAXFORCEPARAM 12
50 #define NR_RBDIHS 6
51 #define NR_FOURDIHS 4
53 typedef atom_id t_iatom;
55 /* this MUST correspond to the
56 t_interaction_function[F_NRE] in gmxlib/ifunc.c */
57 enum {
58 F_BONDS,
59 F_G96BONDS,
60 F_MORSE,
61 F_CUBICBONDS,
62 F_CONNBONDS,
63 F_HARMONIC,
64 F_FENEBONDS,
65 F_TABBONDS,
66 F_TABBONDSNC,
67 F_RESTRBONDS,
68 F_ANGLES,
69 F_G96ANGLES,
70 F_LINEAR_ANGLES,
71 F_CROSS_BOND_BONDS,
72 F_CROSS_BOND_ANGLES,
73 F_UREY_BRADLEY,
74 F_QUARTIC_ANGLES,
75 F_TABANGLES,
76 F_PDIHS,
77 F_RBDIHS,
78 F_FOURDIHS,
79 F_IDIHS,
80 F_PIDIHS,
81 F_TABDIHS,
82 F_CMAP,
83 F_GB12,
84 F_GB13,
85 F_GB14,
86 F_GBPOL,
87 F_NPSOLVATION,
88 F_LJ14,
89 F_COUL14,
90 F_LJC14_Q,
91 F_LJC_PAIRS_NB,
92 F_LJ,
93 F_BHAM,
94 F_LJ_LR,
95 F_BHAM_LR,
96 F_DISPCORR,
97 F_COUL_SR,
98 F_COUL_LR,
99 F_RF_EXCL,
100 F_COUL_RECIP,
101 F_DPD,
102 F_POLARIZATION,
103 F_WATER_POL,
104 F_THOLE_POL,
105 F_ANHARM_POL,
106 F_POSRES,
107 F_DISRES,
108 F_DISRESVIOL,
109 F_ORIRES,
110 F_ORIRESDEV,
111 F_ANGRES,
112 F_ANGRESZ,
113 F_DIHRES,
114 F_DIHRESVIOL,
115 F_CONSTR,
116 F_CONSTRNC,
117 F_SETTLE,
118 F_VSITE2,
119 F_VSITE3,
120 F_VSITE3FD,
121 F_VSITE3FAD,
122 F_VSITE3OUT,
123 F_VSITE4FD,
124 F_VSITE4FDN,
125 F_VSITEN,
126 F_COM_PULL,
127 F_EQM,
128 F_EPOT,
129 F_EKIN,
130 F_ETOT,
131 F_ECONSERVED,
132 F_TEMP,
133 F_VTEMP,
134 F_PDISPCORR,
135 F_PRES,
136 F_DHDL_CON,
137 F_DVDL,
138 F_DKDL,
139 F_DVDL_COUL,
140 F_DVDL_VDW,
141 F_DVDL_BONDED,
142 F_DVDL_RESTRAINT,
143 F_DVDL_TEMPERATURE, /* not calculated for now, but should just be the energy (NVT) or enthalpy (NPT), or 0 (NVE) */
144 F_NRE /* This number is for the total number of energies */
147 #define IS_RESTRAINT_TYPE(ifunc) (((ifunc==F_POSRES) || (ifunc==F_DISRES) || (ifunc==F_RESTRBONDS) || (ifunc==F_DISRESVIOL) || (ifunc==F_ORIRES) || (ifunc==F_ORIRESDEV) || (ifunc==F_ANGRES) || (ifunc == F_ANGRESZ) || (ifunc==F_DIHRES)))
149 /* A macro for checking if ftype is an explicit pair-listed LJ or COULOMB
150 * interaction type:
151 * bonded LJ (usually 1-4), or special listed non-bonded for FEP.
153 #define IS_LISTED_LJ_C(ftype) ((ftype) >= F_LJ14 && (ftype) <= F_LJC_PAIRS_NB)
155 typedef union
157 /* Some parameters have A and B values for free energy calculations.
158 * The B values are not used for regular simulations of course.
159 * Free Energy for nonbondeds can be computed by changing the atom type.
160 * The harmonic type is used for all harmonic potentials:
161 * bonds, angles and improper dihedrals
163 struct {real a,b,c; } bham;
164 struct {real rA,krA,rB,krB; } harmonic;
165 struct {real klinA,aA,klinB,aB; } linangle;
166 struct {real lowA,up1A,up2A,kA,lowB,up1B,up2B,kB; } restraint;
167 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
168 struct {real b0,kb,kcub; } cubic;
169 struct {real bm,kb; } fene;
170 struct {real r1e,r2e,krr; } cross_bb;
171 struct {real r1e,r2e,r3e,krt; } cross_ba;
172 struct {real thetaA,kthetaA,r13A,kUBA,thetaB,kthetaB,r13B,kUBB;} u_b;
173 struct {real theta,c[5]; } qangle;
174 struct {real alpha; } polarize;
175 struct {real alpha,drcut,khyp; } anharm_polarize;
176 struct {real al_x,al_y,al_z,rOH,rHH,rOD; } wpol;
177 struct {real a,alpha1,alpha2,rfac; } thole;
178 struct {real c6,c12; } lj;
179 struct {real c6A,c12A,c6B,c12B; } lj14;
180 struct {real fqq,qi,qj,c6,c12; } ljc14;
181 struct {real qi,qj,c6,c12; } ljcnb;
182 /* Proper dihedrals can not have different multiplicity when
183 * doing free energy calculations, because the potential would not
184 * be periodic anymore.
186 struct {real phiA,cpA;int mult;real phiB,cpB; } pdihs;
187 struct {real dA,dB; } constr;
188 /* Settle can not be used for Free energy calculations of water bond geometry.
189 * Use shake (or lincs) instead if you have to change the water bonds.
191 struct {real doh,dhh; } settle;
192 struct {real b0A,cbA,betaA,b0B,cbB,betaB; } morse;
193 struct {real pos0A[DIM],fcA[DIM],pos0B[DIM],fcB[DIM]; } posres;
194 struct {real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS]; } rbdihs;
195 struct {real a,b,c,d,e,f; } vsite;
196 struct {int n; real a; } vsiten;
197 /* NOTE: npair is only set after reading the tpx file */
198 struct {real low,up1,up2,kfac;int type,label,npair; } disres;
199 struct {real phiA,dphiA,kfacA,phiB,dphiB,kfacB; } dihres;
200 struct {int ex,power,label; real c,obs,kfac; } orires;
201 struct {int table;real kA;real kB; } tab;
202 struct {real sar,st,pi,gbr,bmlt; } gb;
203 struct {int cmapA,cmapB; } cmap;
204 struct {real buf[MAXFORCEPARAM]; } generic; /* Conversion */
205 } t_iparams;
207 typedef int t_functype;
210 * The nonperturbed/perturbed interactions are now separated (sorted) in the
211 * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
212 * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
213 * interactions.
215 typedef struct
217 int nr;
218 int nr_nonperturbed;
219 t_iatom *iatoms;
220 int nalloc;
221 } t_ilist;
224 * The struct t_ilist defines a list of atoms with their interactions.
225 * General field description:
226 * int nr
227 * the size (nr elements) of the interactions array (iatoms[]).
228 * t_iatom *iatoms
229 * specifies which atoms are involved in an interaction of a certain
230 * type. The layout of this array is as follows:
232 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
233 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
234 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
236 * So for interaction type type1 3 atoms are needed, and for type2 and
237 * type3 only 2. The type identifier is used to select the function to
238 * calculate the interaction and its actual parameters. This type
239 * identifier is an index in a params[] and functype[] array.
242 typedef struct
244 real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
245 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
246 } cmapdata_t;
248 typedef struct
250 int ngrid; /* Number of allocated cmap (cmapdata_t ) grids */
251 int grid_spacing; /* Grid spacing */
252 cmapdata_t *cmapdata; /* Pointer to grid with actual, pre-interpolated data */
253 } gmx_cmap_t;
256 typedef struct
258 int ntypes;
259 int atnr;
260 t_functype *functype;
261 t_iparams *iparams;
262 double reppow; /* The repulsion power for VdW: C12*r^-reppow */
263 real fudgeQQ; /* The scaling factor for Coulomb 1-4: f*q1*q2 */
264 gmx_cmap_t cmap_grid; /* The dihedral correction maps */
265 } gmx_ffparams_t;
267 enum {
268 ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
271 typedef struct
273 int ntypes;
274 int atnr;
275 t_functype *functype;
276 t_iparams *iparams;
277 real fudgeQQ;
278 gmx_cmap_t cmap_grid;
279 t_iparams *iparams_posres;
280 int iparams_posres_nalloc;
282 t_ilist il[F_NRE];
283 int ilsort;
284 } t_idef;
287 * The struct t_idef defines all the interactions for the complete
288 * simulation. The structure is setup in such a way that the multinode
289 * version of the program can use it as easy as the single node version.
290 * General field description:
291 * int ntypes
292 * defines the number of elements in functype[] and param[].
293 * int nodeid
294 * the node id (if parallel machines)
295 * int atnr
296 * the number of atomtypes
297 * t_functype *functype
298 * array of length ntypes, defines for every force type what type of
299 * function to use. Every "bond" with the same function but different
300 * force parameters is a different force type. The type identifier in the
301 * forceatoms[] array is an index in this array.
302 * t_iparams *iparams
303 * array of length ntypes, defines the parameters for every interaction
304 * type. The type identifier in the actual interaction list
305 * (ilist[ftype].iatoms[]) is an index in this array.
306 * gmx_cmap_t cmap_grid
307 * the grid for the dihedral pair correction maps.
308 * t_iparams *iparams_posres
309 * defines the parameters for position restraints only.
310 * Position restraints are the only interactions that have different
311 * parameters (reference positions) for different molecules
312 * of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
313 * t_ilist il[F_NRE]
314 * The list of interactions for each type. Note that some,
315 * such as LJ and COUL will have 0 entries.
318 typedef struct {
319 int n; /* n+1 is the number of points */
320 real scale; /* distance between two points */
321 real *tab; /* the actual tables, per point there are 4 numbers */
322 } bondedtable_t;
324 #ifdef __cplusplus
326 #endif
329 #endif