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