Properly finalize MPI on mdrun -version. Fixes #1313
[gromacs.git] / include / types / idef.h
blob5f4ce1e0b7223f75da4fef76f17eebec8690c8e7
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 #ifndef _idef_h
41 #define _idef_h
43 #include "simple.h"
45 #ifdef __cplusplus
46 extern "C" {
47 #endif
50 /* check kernel/toppush.c when you change these numbers */
51 #define MAXATOMLIST 6
52 #define MAXFORCEPARAM 12
53 #define NR_RBDIHS 6
54 #define NR_FOURDIHS 4
56 typedef atom_id t_iatom;
58 /* this MUST correspond to the
59 t_interaction_function[F_NRE] in gmxlib/ifunc.c */
60 enum {
61 F_BONDS,
62 F_G96BONDS,
63 F_MORSE,
64 F_CUBICBONDS,
65 F_CONNBONDS,
66 F_HARMONIC,
67 F_FENEBONDS,
68 F_TABBONDS,
69 F_TABBONDSNC,
70 F_RESTRBONDS,
71 F_ANGLES,
72 F_G96ANGLES,
73 F_LINEAR_ANGLES,
74 F_CROSS_BOND_BONDS,
75 F_CROSS_BOND_ANGLES,
76 F_UREY_BRADLEY,
77 F_QUARTIC_ANGLES,
78 F_TABANGLES,
79 F_PDIHS,
80 F_RBDIHS,
81 F_FOURDIHS,
82 F_IDIHS,
83 F_PIDIHS,
84 F_TABDIHS,
85 F_CMAP,
86 F_GB12,
87 F_GB13,
88 F_GB14,
89 F_GBPOL,
90 F_NPSOLVATION,
91 F_LJ14,
92 F_COUL14,
93 F_LJC14_Q,
94 F_LJC_PAIRS_NB,
95 F_LJ,
96 F_BHAM,
97 F_LJ_LR,
98 F_BHAM_LR,
99 F_DISPCORR,
100 F_COUL_SR,
101 F_COUL_LR,
102 F_RF_EXCL,
103 F_COUL_RECIP,
104 F_DPD,
105 F_POLARIZATION,
106 F_WATER_POL,
107 F_THOLE_POL,
108 F_ANHARM_POL,
109 F_POSRES,
110 F_DISRES,
111 F_DISRESVIOL,
112 F_ORIRES,
113 F_ORIRESDEV,
114 F_ANGRES,
115 F_ANGRESZ,
116 F_DIHRES,
117 F_DIHRESVIOL,
118 F_CONSTR,
119 F_CONSTRNC,
120 F_SETTLE,
121 F_VSITE2,
122 F_VSITE3,
123 F_VSITE3FD,
124 F_VSITE3FAD,
125 F_VSITE3OUT,
126 F_VSITE4FD,
127 F_VSITE4FDN,
128 F_VSITEN,
129 F_COM_PULL,
130 F_EQM,
131 F_EPOT,
132 F_EKIN,
133 F_ETOT,
134 F_ECONSERVED,
135 F_TEMP,
136 F_VTEMP_NOLONGERUSED,
137 F_PDISPCORR,
138 F_PRES,
139 F_DVDL_CONSTR,
140 F_DVDL,
141 F_DKDL,
142 F_DVDL_COUL,
143 F_DVDL_VDW,
144 F_DVDL_BONDED,
145 F_DVDL_RESTRAINT,
146 F_DVDL_TEMPERATURE, /* not calculated for now, but should just be the energy (NVT) or enthalpy (NPT), or 0 (NVE) */
147 F_NRE /* This number is for the total number of energies */
150 #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)))
152 /* A macro for checking if ftype is an explicit pair-listed LJ or COULOMB
153 * interaction type:
154 * bonded LJ (usually 1-4), or special listed non-bonded for FEP.
156 #define IS_LISTED_LJ_C(ftype) ((ftype) >= F_LJ14 && (ftype) <= F_LJC_PAIRS_NB)
158 typedef union
160 /* Some parameters have A and B values for free energy calculations.
161 * The B values are not used for regular simulations of course.
162 * Free Energy for nonbondeds can be computed by changing the atom type.
163 * The harmonic type is used for all harmonic potentials:
164 * bonds, angles and improper dihedrals
166 struct {
167 real a, b, c;
168 } bham;
169 struct {
170 real rA, krA, rB, krB;
171 } harmonic;
172 struct {
173 real klinA, aA, klinB, aB;
174 } linangle;
175 struct {
176 real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
177 } restraint;
178 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
179 struct {
180 real b0, kb, kcub;
181 } cubic;
182 struct {
183 real bm, kb;
184 } fene;
185 struct {
186 real r1e, r2e, krr;
187 } cross_bb;
188 struct {
189 real r1e, r2e, r3e, krt;
190 } cross_ba;
191 struct {
192 real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
193 } u_b;
194 struct {
195 real theta, c[5];
196 } qangle;
197 struct {
198 real alpha;
199 } polarize;
200 struct {
201 real alpha, drcut, khyp;
202 } anharm_polarize;
203 struct {
204 real al_x, al_y, al_z, rOH, rHH, rOD;
205 } wpol;
206 struct {
207 real a, alpha1, alpha2, rfac;
208 } thole;
209 struct {
210 real c6, c12;
211 } lj;
212 struct {
213 real c6A, c12A, c6B, c12B;
214 } lj14;
215 struct {
216 real fqq, qi, qj, c6, c12;
217 } ljc14;
218 struct {
219 real qi, qj, c6, c12;
220 } ljcnb;
221 /* Proper dihedrals can not have different multiplicity when
222 * doing free energy calculations, because the potential would not
223 * be periodic anymore.
225 struct {
226 real phiA, cpA; int mult; real phiB, cpB;
227 } pdihs;
228 struct {
229 real dA, dB;
230 } constr;
231 /* Settle can not be used for Free energy calculations of water bond geometry.
232 * Use shake (or lincs) instead if you have to change the water bonds.
234 struct {
235 real doh, dhh;
236 } settle;
237 struct {
238 real b0A, cbA, betaA, b0B, cbB, betaB;
239 } morse;
240 struct {
241 real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
242 } posres;
243 struct {
244 real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
245 } rbdihs;
246 struct {
247 real a, b, c, d, e, f;
248 } vsite;
249 struct {
250 int n; real a;
251 } vsiten;
252 /* NOTE: npair is only set after reading the tpx file */
253 struct {
254 real low, up1, up2, kfac; int type, label, npair;
255 } disres;
256 struct {
257 real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
258 } dihres;
259 struct {
260 int ex, power, label; real c, obs, kfac;
261 } orires;
262 struct {
263 int table; real kA; real kB;
264 } tab;
265 struct {
266 real sar, st, pi, gbr, bmlt;
267 } gb;
268 struct {
269 int cmapA, cmapB;
270 } cmap;
271 struct {
272 real buf[MAXFORCEPARAM];
273 } generic; /* Conversion */
274 } t_iparams;
276 typedef int t_functype;
279 * The nonperturbed/perturbed interactions are now separated (sorted) in the
280 * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
281 * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
282 * interactions.
284 typedef struct
286 int nr;
287 int nr_nonperturbed;
288 t_iatom *iatoms;
289 int nalloc;
290 } t_ilist;
293 * The struct t_ilist defines a list of atoms with their interactions.
294 * General field description:
295 * int nr
296 * the size (nr elements) of the interactions array (iatoms[]).
297 * t_iatom *iatoms
298 * specifies which atoms are involved in an interaction of a certain
299 * type. The layout of this array is as follows:
301 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
302 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
303 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
305 * So for interaction type type1 3 atoms are needed, and for type2 and
306 * type3 only 2. The type identifier is used to select the function to
307 * calculate the interaction and its actual parameters. This type
308 * identifier is an index in a params[] and functype[] array.
311 typedef struct
313 real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
314 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
315 } cmapdata_t;
317 typedef struct
319 int ngrid; /* Number of allocated cmap (cmapdata_t ) grids */
320 int grid_spacing; /* Grid spacing */
321 cmapdata_t *cmapdata; /* Pointer to grid with actual, pre-interpolated data */
322 } gmx_cmap_t;
325 typedef struct
327 int ntypes;
328 int atnr;
329 t_functype *functype;
330 t_iparams *iparams;
331 double reppow; /* The repulsion power for VdW: C12*r^-reppow */
332 real fudgeQQ; /* The scaling factor for Coulomb 1-4: f*q1*q2 */
333 gmx_cmap_t cmap_grid; /* The dihedral correction maps */
334 } gmx_ffparams_t;
336 enum {
337 ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
340 typedef struct
342 int ntypes;
343 int atnr;
344 t_functype *functype;
345 t_iparams *iparams;
346 real fudgeQQ;
347 gmx_cmap_t cmap_grid;
348 t_iparams *iparams_posres;
349 int iparams_posres_nalloc;
351 t_ilist il[F_NRE];
352 int ilsort;
353 } t_idef;
356 * The struct t_idef defines all the interactions for the complete
357 * simulation. The structure is setup in such a way that the multinode
358 * version of the program can use it as easy as the single node version.
359 * General field description:
360 * int ntypes
361 * defines the number of elements in functype[] and param[].
362 * int nodeid
363 * the node id (if parallel machines)
364 * int atnr
365 * the number of atomtypes
366 * t_functype *functype
367 * array of length ntypes, defines for every force type what type of
368 * function to use. Every "bond" with the same function but different
369 * force parameters is a different force type. The type identifier in the
370 * forceatoms[] array is an index in this array.
371 * t_iparams *iparams
372 * array of length ntypes, defines the parameters for every interaction
373 * type. The type identifier in the actual interaction list
374 * (ilist[ftype].iatoms[]) is an index in this array.
375 * gmx_cmap_t cmap_grid
376 * the grid for the dihedral pair correction maps.
377 * t_iparams *iparams_posres
378 * defines the parameters for position restraints only.
379 * Position restraints are the only interactions that have different
380 * parameters (reference positions) for different molecules
381 * of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
382 * t_ilist il[F_NRE]
383 * The list of interactions for each type. Note that some,
384 * such as LJ and COUL will have 0 entries.
387 typedef struct {
388 int n; /* n+1 is the number of points */
389 real scale; /* distance between two points */
390 real *data; /* the actual table data, per point there are 4 numbers */
391 } bondedtable_t;
393 #ifdef __cplusplus
395 #endif
398 #endif