Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / include / types / idef.h
blobefa74f8059fd1c06e38955204ea18c0b28cccfe0
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_CROSS_BOND_BONDS,
71 F_CROSS_BOND_ANGLES,
72 F_UREY_BRADLEY,
73 F_QUARTIC_ANGLES,
74 F_TABANGLES,
75 F_PDIHS,
76 F_RBDIHS,
77 F_FOURDIHS,
78 F_IDIHS,
79 F_PIDIHS,
80 F_TABDIHS,
81 F_CMAP,
82 F_GB12,
83 F_GB13,
84 F_GB14,
85 F_GBPOL,
86 F_NPSOLVATION,
87 F_LJ14,
88 F_COUL14,
89 F_LJC14_Q,
90 F_LJC_PAIRS_NB,
91 F_LJ,
92 F_BHAM,
93 F_LJ_LR,
94 F_BHAM_LR,
95 F_DISPCORR,
96 F_COUL_SR,
97 F_COUL_LR,
98 F_RF_EXCL,
99 F_COUL_RECIP,
100 F_DPD,
101 F_POLARIZATION,
102 F_WATER_POL,
103 F_THOLE_POL,
104 F_POSRES,
105 F_DISRES,
106 F_DISRESVIOL,
107 F_ORIRES,
108 F_ORIRESDEV,
109 F_ANGRES,
110 F_ANGRESZ,
111 F_DIHRES,
112 F_DIHRESVIOL,
113 F_CONSTR,
114 F_CONSTRNC,
115 F_SETTLE,
116 F_VSITE2,
117 F_VSITE3,
118 F_VSITE3FD,
119 F_VSITE3FAD,
120 F_VSITE3OUT,
121 F_VSITE4FD,
122 F_VSITE4FDN,
123 F_VSITEN,
124 F_COM_PULL,
125 F_EQM,
126 F_EPOT,
127 F_EKIN,
128 F_ETOT,
129 F_ECONSERVED,
130 F_TEMP,
131 F_VTEMP,
132 F_PDISPCORR,
133 F_PRES,
134 F_DVDL,
135 F_DKDL,
136 F_DHDL_CON,
137 F_NRE /* This number is for the total number of energies */
140 typedef union
142 /* Some parameters have A and B values for free energy calculations.
143 * The B values are not used for regular simulations of course.
144 * Free Energy for nonbondeds can be computed by changing the atom type.
145 * The harmonic type is used for all harmonic potentials:
146 * bonds, angles and improper dihedrals
148 struct {real a,b,c; } bham;
149 struct {real rA,krA,rB,krB; } harmonic;
150 struct {real lowA,up1A,up2A,kA,lowB,up1B,up2B,kB; } restraint;
151 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
152 struct {real b0,kb,kcub; } cubic;
153 struct {real bm,kb; } fene;
154 struct {real r1e,r2e,krr; } cross_bb;
155 struct {real r1e,r2e,r3e,krt; } cross_ba;
156 struct {real theta,ktheta,r13,kUB; } u_b;
157 struct {real theta,c[5]; } qangle;
158 struct {real alpha; } polarize;
159 struct {real al_x,al_y,al_z,rOH,rHH,rOD; } wpol;
160 struct {real a,alpha1,alpha2,rfac; } thole;
161 struct {real c6,c12; } lj;
162 struct {real c6A,c12A,c6B,c12B; } lj14;
163 struct {real fqq,qi,qj,c6,c12; } ljc14;
164 struct {real qi,qj,c6,c12; } ljcnb;
165 /* Proper dihedrals can not have different multiplicity when
166 * doing free energy calculations, because the potential would not
167 * be periodic anymore.
169 struct {real phiA,cpA;int mult;real phiB,cpB; } pdihs;
170 struct {real dA,dB; } constr;
171 /* Settle can not be used for Free energy calculations of water bond geometry.
172 * Use shake (or lincs) instead if you have to change the water bonds.
174 struct {real doh,dhh; } settle;
175 /* No free energy supported for morse bonds */
176 struct {real b0,cb,beta; } morse;
177 struct {real pos0A[DIM],fcA[DIM],pos0B[DIM],fcB[DIM]; } posres;
178 struct {real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS]; } rbdihs;
179 struct {real a,b,c,d,e,f; } vsite;
180 struct {int n; real a; } vsiten;
181 /* NOTE: npair is only set after reading the tpx file */
182 struct {real low,up1,up2,kfac;int type,label,npair; } disres;
183 struct {real phi,dphi,kfac;int label,power; } dihres;
184 struct {int ex,power,label; real c,obs,kfac; } orires;
185 struct {int table;real kA;real kB; } tab;
186 struct {real sar,st,pi,gbr,bmlt; } gb;
187 struct {int cmapA,cmapB; } cmap;
188 struct {real buf[MAXFORCEPARAM]; } generic; /* Conversion */
189 } t_iparams;
191 typedef int t_functype;
194 * The nonperturbed/perturbed interactions are now separated (sorted) in the
195 * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
196 * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
197 * interactions.
199 typedef struct
201 int nr;
202 int nr_nonperturbed;
203 t_iatom *iatoms;
204 int nalloc;
205 } t_ilist;
208 * The struct t_ilist defines a list of atoms with their interactions.
209 * General field description:
210 * int nr
211 * the size (nr elements) of the interactions array (iatoms[]).
212 * t_iatom *iatoms
213 * specifies which atoms are involved in an interaction of a certain
214 * type. The layout of this array is as follows:
216 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
217 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
218 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
220 * So for interaction type type1 3 atoms are needed, and for type2 and
221 * type3 only 2. The type identifier is used to select the function to
222 * calculate the interaction and its actual parameters. This type
223 * identifier is an index in a params[] and functype[] array.
226 typedef struct
228 real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
229 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
230 } cmapdata_t;
232 typedef struct
234 int ngrid; /* Number of allocated cmap (cmapdata_t ) grids */
235 int grid_spacing; /* Grid spacing */
236 cmapdata_t *cmapdata; /* Pointer to grid with actual, pre-interpolated data */
237 } gmx_cmap_t;
240 typedef struct
242 int ntypes;
243 int atnr;
244 t_functype *functype;
245 t_iparams *iparams;
246 double reppow; /* The repulsion power for VdW: C12*r^-reppow */
247 real fudgeQQ; /* The scaling factor for Coulomb 1-4: f*q1*q2 */
248 gmx_cmap_t cmap_grid; /* The dihedral correction maps */
249 } gmx_ffparams_t;
251 enum {
252 ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
255 typedef struct
257 int ntypes;
258 int atnr;
259 t_functype *functype;
260 t_iparams *iparams;
261 real fudgeQQ;
262 gmx_cmap_t cmap_grid;
263 t_iparams *iparams_posres;
264 int iparams_posres_nalloc;
266 t_ilist il[F_NRE];
267 int ilsort;
268 } t_idef;
271 * The struct t_idef defines all the interactions for the complete
272 * simulation. The structure is setup in such a way that the multinode
273 * version of the program can use it as easy as the single node version.
274 * General field description:
275 * int ntypes
276 * defines the number of elements in functype[] and param[].
277 * int nodeid
278 * the node id (if parallel machines)
279 * int atnr
280 * the number of atomtypes
281 * t_functype *functype
282 * array of length ntypes, defines for every force type what type of
283 * function to use. Every "bond" with the same function but different
284 * force parameters is a different force type. The type identifier in the
285 * forceatoms[] array is an index in this array.
286 * t_iparams *iparams
287 * array of length ntypes, defines the parameters for every interaction
288 * type. The type identifier in the actual interaction list
289 * (ilist[ftype].iatoms[]) is an index in this array.
290 * gmx_cmap_t cmap_grid
291 * the grid for the dihedral pair correction maps.
292 * t_iparams *iparams_posres
293 * defines the parameters for position restraints only.
294 * Position restraints are the only interactions that have different
295 * parameters (reference positions) for different molecules
296 * of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
297 * t_ilist il[F_NRE]
298 * The list of interactions for each type. Note that some,
299 * such as LJ and COUL will have 0 entries.
302 typedef struct {
303 int n; /* n+1 is the number of points */
304 real scale; /* distance between two points */
305 real *tab; /* the actual tables, per point there are 4 numbers */
306 } bondedtable_t;
308 #ifdef __cplusplus
310 #endif
313 #endif