Move pair-interaction code into listed-forces module
[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 typedef union t_iparams
156 {
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
162      */
163     struct {
164         real a, b, c;
165     } bham;
166     struct {
167         real rA, krA, rB, krB;
168     } harmonic;
169     struct {
170         real klinA, aA, klinB, aB;
171     } linangle;
172     struct {
173         real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
174     } restraint;
175     /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
176     struct {
177         real b0, kb, kcub;
178     } cubic;
179     struct {
180         real bm, kb;
181     } fene;
182     struct {
183         real r1e, r2e, krr;
184     } cross_bb;
185     struct {
186         real r1e, r2e, r3e, krt;
187     } cross_ba;
188     struct {
189         real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
190     } u_b;
191     struct {
192         real theta, c[5];
193     } qangle;
194     struct {
195         real alpha;
196     } polarize;
197     struct {
198         real alpha, drcut, khyp;
199     } anharm_polarize;
200     struct {
201         real al_x, al_y, al_z, rOH, rHH, rOD;
202     } wpol;
203     struct {
204         real a, alpha1, alpha2, rfac;
205     } thole;
206     struct {
207         real c6, c12;
208     } lj;
209     struct {
210         real c6A, c12A, c6B, c12B;
211     } lj14;
212     struct {
213         real fqq, qi, qj, c6, c12;
214     } ljc14;
215     struct {
216         real qi, qj, c6, c12;
217     } ljcnb;
218     /* Proper dihedrals can not have different multiplicity when
219      * doing free energy calculations, because the potential would not
220      * be periodic anymore.
221      */
222     struct {
223         real phiA, cpA; int mult; real phiB, cpB;
224     } pdihs;
225     struct {
226         real dA, dB;
227     } constr;
228     /* Settle can not be used for Free energy calculations of water bond geometry.
229      * Use shake (or lincs) instead if you have to change the water bonds.
230      */
231     struct {
232         real doh, dhh;
233     } settle;
234     struct {
235         real b0A, cbA, betaA, b0B, cbB, betaB;
236     } morse;
237     struct {
238         real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
239     } posres;
240     struct {
241         real pos0[DIM], r, k; int geom;
242     } fbposres;
243     struct {
244         real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
245     } rbdihs;
246     struct {
247         real cbtcA[NR_CBTDIHS], cbtcB[NR_CBTDIHS];
248     } cbtdihs;
249     struct {
250         real a, b, c, d, e, f;
251     } vsite;
252     struct {
253         int  n; real a;
254     } vsiten;
255     /* NOTE: npair is only set after reading the tpx file */
256     struct {
257         real low, up1, up2, kfac; int type, label, npair;
258     } disres;
259     struct {
260         real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
261     } dihres;
262     struct {
263         int  ex, power, label; real c, obs, kfac;
264     } orires;
265     struct {
266         int  table; real kA; real kB;
267     } tab;
268     struct {
269         real sar, st, pi, gbr, bmlt;
270     } gb;
271     struct {
272         int cmapA, cmapB;
273     } cmap;
274     struct {
275         real buf[MAXFORCEPARAM];
276     } generic;                                               /* Conversion */
277 } t_iparams;
278
279 typedef int t_functype;
280
281 /*
282  * The nonperturbed/perturbed interactions are now separated (sorted) in the
283  * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
284  * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
285  * interactions.
286  */
287 typedef struct t_ilist
288 {
289     int      nr;
290     int      nr_nonperturbed;
291     t_iatom *iatoms;
292     int      nalloc;
293 } t_ilist;
294
295 /*
296  * The struct t_ilist defines a list of atoms with their interactions.
297  * General field description:
298  *   int nr
299  *      the size (nr elements) of the interactions array (iatoms[]).
300  *   t_iatom *iatoms
301  *  specifies which atoms are involved in an interaction of a certain
302  *       type. The layout of this array is as follows:
303  *
304  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
305  *        |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
306  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
307  *
308  *  So for interaction type type1 3 atoms are needed, and for type2 and
309  *      type3 only 2. The type identifier is used to select the function to
310  *      calculate the interaction and its actual parameters. This type
311  *      identifier is an index in a params[] and functype[] array.
312  */
313
314 typedef struct
315 {
316     real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
317     /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
318 } gmx_cmapdata_t;
319
320 typedef struct gmx_cmap_t
321 {
322     int             ngrid;        /* Number of allocated cmap (cmapdata_t ) grids */
323     int             grid_spacing; /* Grid spacing */
324     gmx_cmapdata_t *cmapdata;     /* Pointer to grid with actual, pre-interpolated data */
325 } gmx_cmap_t;
326
327
328 typedef struct gmx_ffparams_t
329 {
330     int         ntypes;
331     int         atnr;
332     t_functype *functype;
333     t_iparams  *iparams;
334     double      reppow;    /* The repulsion power for VdW: C12*r^-reppow   */
335     real        fudgeQQ;   /* The scaling factor for Coulomb 1-4: f*q1*q2  */
336     gmx_cmap_t  cmap_grid; /* The dihedral correction maps                 */
337 } gmx_ffparams_t;
338
339 enum {
340     ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
341 };
342
343 typedef struct t_idef
344 {
345     int         ntypes;
346     int         atnr;
347     t_functype *functype;
348     t_iparams  *iparams;
349     real        fudgeQQ;
350     gmx_cmap_t  cmap_grid;
351     t_iparams  *iparams_posres, *iparams_fbposres;
352     int         iparams_posres_nalloc, iparams_fbposres_nalloc;
353
354     t_ilist     il[F_NRE];
355     int         ilsort;
356     int         nthreads;
357     int        *il_thread_division;
358     int         il_thread_division_nalloc;
359 } t_idef;
360
361 /*
362  * The struct t_idef defines all the interactions for the complete
363  * simulation. The structure is setup in such a way that the multinode
364  * version of the program  can use it as easy as the single node version.
365  * General field description:
366  *   int ntypes
367  *      defines the number of elements in functype[] and param[].
368  *   int nodeid
369  *      the node id (if parallel machines)
370  *   int atnr
371  *      the number of atomtypes
372  *   t_functype *functype
373  *      array of length ntypes, defines for every force type what type of
374  *      function to use. Every "bond" with the same function but different
375  *      force parameters is a different force type. The type identifier in the
376  *      forceatoms[] array is an index in this array.
377  *   t_iparams *iparams
378  *      array of length ntypes, defines the parameters for every interaction
379  *      type. The type identifier in the actual interaction list
380  *      (ilist[ftype].iatoms[]) is an index in this array.
381  *   gmx_cmap_t cmap_grid
382  *      the grid for the dihedral pair correction maps.
383  *   t_iparams *iparams_posres, *iparams_fbposres
384  *      defines the parameters for position restraints only.
385  *      Position restraints are the only interactions that have different
386  *      parameters (reference positions) for different molecules
387  *      of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
388  *   t_ilist il[F_NRE]
389  *      The list of interactions for each type. Note that some,
390  *      such as LJ and COUL will have 0 entries.
391  *   int ilsort
392  *      The state of the sorting of il, values are provided above.
393  *   int nthreads
394  *      The number of threads used to set il_thread_division.
395  *   int *il_thread_division
396  *      The division of the normal bonded interactions of threads.
397  *      il_thread_division[ftype*(nthreads+1)+t] contains an index
398  *      into il[ftype].iatoms; thread th operates on t=th to t=th+1.
399  *   int il_thread_division_nalloc
400  *      The allocated size of il_thread_division,
401  *      should be at least F_NRE*(nthreads+1).
402  */
403
404 #ifdef __cplusplus
405 }
406 #endif
407
408 #endif