Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / legacyheaders / types / idef.h
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GRoups of Organic Molecules in ACtion for Science
34  */
35
36
37 #ifndef _idef_h
38 #define _idef_h
39
40 #include "simple.h"
41
42 #ifdef __cplusplus
43 extern "C" {
44 #endif
45
46
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
52
53 typedef atom_id t_iatom;
54
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_FBPOSRES,
108     F_DISRES,
109     F_DISRESVIOL,
110     F_ORIRES,
111     F_ORIRESDEV,
112     F_ANGRES,
113     F_ANGRESZ,
114     F_DIHRES,
115     F_DIHRESVIOL,
116     F_CONSTR,
117     F_CONSTRNC,
118     F_SETTLE,
119     F_VSITE2,
120     F_VSITE3,
121     F_VSITE3FD,
122     F_VSITE3FAD,
123     F_VSITE3OUT,
124     F_VSITE4FD,
125     F_VSITE4FDN,
126     F_VSITEN,
127     F_COM_PULL,
128     F_EQM,
129     F_EPOT,
130     F_EKIN,
131     F_ETOT,
132     F_ECONSERVED,
133     F_TEMP,
134     F_VTEMP_NOLONGERUSED,
135     F_PDISPCORR,
136     F_PRES,
137     F_DVDL_CONSTR,
138     F_DVDL,
139     F_DKDL,
140     F_DVDL_COUL,
141     F_DVDL_VDW,
142     F_DVDL_BONDED,
143     F_DVDL_RESTRAINT,
144     F_DVDL_TEMPERATURE, /* not calculated for now, but should just be the energy (NVT) or enthalpy (NPT), or 0 (NVE) */
145     F_NRE               /* This number is for the total number of energies      */
146 };
147
148 #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
150 /* A macro for checking if ftype is an explicit pair-listed LJ or COULOMB
151  * interaction type:
152  * bonded LJ (usually 1-4), or special listed non-bonded for FEP.
153  */
154 #define IS_LISTED_LJ_C(ftype) ((ftype) >= F_LJ14 && (ftype) <= F_LJC_PAIRS_NB)
155
156 typedef union
157 {
158     /* Some parameters have A and B values for free energy calculations.
159      * The B values are not used for regular simulations of course.
160      * Free Energy for nonbondeds can be computed by changing the atom type.
161      * The harmonic type is used for all harmonic potentials:
162      * bonds, angles and improper dihedrals
163      */
164     struct {
165         real a, b, c;
166     } bham;
167     struct {
168         real rA, krA, rB, krB;
169     } harmonic;
170     struct {
171         real klinA, aA, klinB, aB;
172     } linangle;
173     struct {
174         real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
175     } restraint;
176     /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
177     struct {
178         real b0, kb, kcub;
179     } cubic;
180     struct {
181         real bm, kb;
182     } fene;
183     struct {
184         real r1e, r2e, krr;
185     } cross_bb;
186     struct {
187         real r1e, r2e, r3e, krt;
188     } cross_ba;
189     struct {
190         real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
191     } u_b;
192     struct {
193         real theta, c[5];
194     } qangle;
195     struct {
196         real alpha;
197     } polarize;
198     struct {
199         real alpha, drcut, khyp;
200     } anharm_polarize;
201     struct {
202         real al_x, al_y, al_z, rOH, rHH, rOD;
203     } wpol;
204     struct {
205         real a, alpha1, alpha2, rfac;
206     } thole;
207     struct {
208         real c6, c12;
209     } lj;
210     struct {
211         real c6A, c12A, c6B, c12B;
212     } lj14;
213     struct {
214         real fqq, qi, qj, c6, c12;
215     } ljc14;
216     struct {
217         real qi, qj, c6, c12;
218     } ljcnb;
219     /* Proper dihedrals can not have different multiplicity when
220      * doing free energy calculations, because the potential would not
221      * be periodic anymore.
222      */
223     struct {
224         real phiA, cpA; int mult; real phiB, cpB;
225     } pdihs;
226     struct {
227         real dA, dB;
228     } constr;
229     /* Settle can not be used for Free energy calculations of water bond geometry.
230      * Use shake (or lincs) instead if you have to change the water bonds.
231      */
232     struct {
233         real doh, dhh;
234     } settle;
235     struct {
236         real b0A, cbA, betaA, b0B, cbB, betaB;
237     } morse;
238     struct {
239         real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
240     } posres;
241     struct {
242         real pos0[DIM], r, k; int geom;
243     } fbposres;
244     struct {
245         real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
246     } rbdihs;
247     struct {
248         real a, b, c, d, e, f;
249     } vsite;
250     struct {
251         int  n; real a;
252     } vsiten;
253     /* NOTE: npair is only set after reading the tpx file */
254     struct {
255         real low, up1, up2, kfac; int type, label, npair;
256     } disres;
257     struct {
258         real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
259     } dihres;
260     struct {
261         int  ex, power, label; real c, obs, kfac;
262     } orires;
263     struct {
264         int  table; real kA; real kB;
265     } tab;
266     struct {
267         real sar, st, pi, gbr, bmlt;
268     } gb;
269     struct {
270         int cmapA, cmapB;
271     } cmap;
272     struct {
273         real buf[MAXFORCEPARAM];
274     } generic;                                               /* Conversion */
275 } t_iparams;
276
277 typedef int t_functype;
278
279 /*
280  * The nonperturbed/perturbed interactions are now separated (sorted) in the
281  * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
282  * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
283  * interactions.
284  */
285 typedef struct
286 {
287     int      nr;
288     int      nr_nonperturbed;
289     t_iatom *iatoms;
290     int      nalloc;
291 } t_ilist;
292
293 /*
294  * The struct t_ilist defines a list of atoms with their interactions.
295  * General field description:
296  *   int nr
297  *      the size (nr elements) of the interactions array (iatoms[]).
298  *   t_iatom *iatoms
299  *  specifies which atoms are involved in an interaction of a certain
300  *       type. The layout of this array is as follows:
301  *
302  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
303  *        |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
304  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
305  *
306  *  So for interaction type type1 3 atoms are needed, and for type2 and
307  *      type3 only 2. The type identifier is used to select the function to
308  *      calculate the interaction and its actual parameters. This type
309  *      identifier is an index in a params[] and functype[] array.
310  */
311
312 typedef struct
313 {
314     real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
315     /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
316 } cmapdata_t;
317
318 typedef struct
319 {
320     int         ngrid;        /* Number of allocated cmap (cmapdata_t ) grids */
321     int         grid_spacing; /* Grid spacing */
322     cmapdata_t *cmapdata;     /* Pointer to grid with actual, pre-interpolated data */
323 } gmx_cmap_t;
324
325
326 typedef struct
327 {
328     int         ntypes;
329     int         atnr;
330     t_functype *functype;
331     t_iparams  *iparams;
332     double      reppow;    /* The repulsion power for VdW: C12*r^-reppow   */
333     real        fudgeQQ;   /* The scaling factor for Coulomb 1-4: f*q1*q2  */
334     gmx_cmap_t  cmap_grid; /* The dihedral correction maps                 */
335 } gmx_ffparams_t;
336
337 enum {
338     ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
339 };
340
341 typedef struct
342 {
343     int         ntypes;
344     int         atnr;
345     t_functype *functype;
346     t_iparams  *iparams;
347     real        fudgeQQ;
348     gmx_cmap_t  cmap_grid;
349     t_iparams  *iparams_posres, *iparams_fbposres;
350     int         iparams_posres_nalloc, iparams_fbposres_nalloc;
351
352     t_ilist     il[F_NRE];
353     int         ilsort;
354 } t_idef;
355
356 /*
357  * The struct t_idef defines all the interactions for the complete
358  * simulation. The structure is setup in such a way that the multinode
359  * version of the program  can use it as easy as the single node version.
360  * General field description:
361  *   int ntypes
362  *      defines the number of elements in functype[] and param[].
363  *   int nodeid
364  *      the node id (if parallel machines)
365  *   int atnr
366  *      the number of atomtypes
367  *   t_functype *functype
368  *      array of length ntypes, defines for every force type what type of
369  *      function to use. Every "bond" with the same function but different
370  *      force parameters is a different force type. The type identifier in the
371  *      forceatoms[] array is an index in this array.
372  *   t_iparams *iparams
373  *      array of length ntypes, defines the parameters for every interaction
374  *      type. The type identifier in the actual interaction list
375  *      (ilist[ftype].iatoms[]) is an index in this array.
376  *   gmx_cmap_t cmap_grid
377  *      the grid for the dihedral pair correction maps.
378  *   t_iparams *iparams_posres, *iparams_fbposres
379  *      defines the parameters for position restraints only.
380  *      Position restraints are the only interactions that have different
381  *      parameters (reference positions) for different molecules
382  *      of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
383  *   t_ilist il[F_NRE]
384  *      The list of interactions for each type. Note that some,
385  *      such as LJ and COUL will have 0 entries.
386  */
387
388 typedef struct {
389     int   n;      /* n+1 is the number of points */
390     real  scale;  /* distance between two points */
391     real *data;   /* the actual table data, per point there are 4 numbers */
392 } bondedtable_t;
393
394 #ifdef __cplusplus
395 }
396 #endif
397
398
399 #endif