Fixing copyright issues and code contributors
[alexxy/gromacs.git] / include / types / 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  * 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38
39
40 #ifndef _idef_h
41 #define _idef_h
42
43 #include "simple.h"
44
45 #ifdef __cplusplus
46 extern "C" {
47 #endif
48
49
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
55
56 typedef atom_id t_iatom;
57
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_DHDL_CON,
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      */
148 };
149
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)))
151
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.
155  */
156 #define IS_LISTED_LJ_C(ftype) ((ftype) >= F_LJ14 && (ftype) <= F_LJC_PAIRS_NB)
157
158 typedef union
159 {
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
165    */
166   struct {real a,b,c;                                      } bham;
167   struct {real rA,krA,rB,krB;                              } harmonic;
168   struct {real klinA,aA,klinB,aB;                          } linangle;
169   struct {real lowA,up1A,up2A,kA,lowB,up1B,up2B,kB;        } restraint;
170   /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */ 
171   struct {real b0,kb,kcub;                                 } cubic;
172   struct {real bm,kb;                                      } fene;
173   struct {real r1e,r2e,krr;                                } cross_bb;
174   struct {real r1e,r2e,r3e,krt;                            } cross_ba;
175   struct {real thetaA,kthetaA,r13A,kUBA,thetaB,kthetaB,r13B,kUBB;} u_b;
176   struct {real theta,c[5];                                 } qangle; 
177   struct {real alpha;                                      } polarize;
178   struct {real alpha,drcut,khyp;                           } anharm_polarize;
179   struct {real al_x,al_y,al_z,rOH,rHH,rOD;                 } wpol;
180   struct {real a,alpha1,alpha2,rfac;                       } thole;
181   struct {real c6,c12;                                     } lj;
182   struct {real c6A,c12A,c6B,c12B;                          } lj14;
183   struct {real fqq,qi,qj,c6,c12;                           } ljc14;
184   struct {real qi,qj,c6,c12;                               } ljcnb;
185   /* Proper dihedrals can not have different multiplicity when
186    * doing free energy calculations, because the potential would not
187    * be periodic anymore.
188    */ 
189   struct {real phiA,cpA;int mult;real phiB,cpB;            } pdihs;
190   struct {real dA,dB;                                      } constr;
191   /* Settle can not be used for Free energy calculations of water bond geometry.
192    * Use shake (or lincs) instead if you have to change the water bonds.
193    */
194   struct {real doh,dhh;                                   } settle;
195   struct {real b0A,cbA,betaA,b0B,cbB,betaB;               } morse;
196   struct {real pos0A[DIM],fcA[DIM],pos0B[DIM],fcB[DIM];   } posres;
197   struct {real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];          } rbdihs;
198   struct {real a,b,c,d,e,f;                               } vsite;   
199   struct {int  n; real a;                                 } vsiten;   
200   /* NOTE: npair is only set after reading the tpx file */
201   struct {real low,up1,up2,kfac;int type,label,npair;     } disres; 
202   struct {real phiA,dphiA,kfacA,phiB,dphiB,kfacB;         } dihres;
203   struct {int  ex,power,label; real c,obs,kfac;           } orires;
204   struct {int  table;real kA;real kB;                     } tab;
205   struct {real sar,st,pi,gbr,bmlt;                        } gb;
206   struct {int cmapA,cmapB;                                } cmap;
207   struct {real buf[MAXFORCEPARAM];                        } generic; /* Conversion */
208 } t_iparams;
209
210 typedef int t_functype;
211
212 /*
213  * The nonperturbed/perturbed interactions are now separated (sorted) in the
214  * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and 
215  * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded 
216  * interactions.
217  */
218 typedef struct
219 {
220   int nr;
221   int nr_nonperturbed;
222   t_iatom *iatoms;
223   int nalloc;
224 } t_ilist;
225
226 /*
227  * The struct t_ilist defines a list of atoms with their interactions. 
228  * General field description:
229  *   int nr
230  *      the size (nr elements) of the interactions array (iatoms[]).
231  *   t_iatom *iatoms
232  *      specifies which atoms are involved in an interaction of a certain 
233  *       type. The layout of this array is as follows:
234  *
235  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
236  *        |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
237  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
238  *
239  *      So for interaction type type1 3 atoms are needed, and for type2 and 
240  *      type3 only 2. The type identifier is used to select the function to 
241  *      calculate the interaction and its actual parameters. This type 
242  *      identifier is an index in a params[] and functype[] array.
243  */
244
245 typedef struct
246 {
247         real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
248         /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
249 } cmapdata_t;
250
251 typedef struct
252 {
253         int ngrid;            /* Number of allocated cmap (cmapdata_t ) grids */
254         int grid_spacing;     /* Grid spacing */
255         cmapdata_t *cmapdata; /* Pointer to grid with actual, pre-interpolated data */
256 } gmx_cmap_t;
257
258
259 typedef struct
260 {
261   int        ntypes;
262   int        atnr;
263   t_functype *functype;
264   t_iparams  *iparams;
265   double     reppow;     /* The repulsion power for VdW: C12*r^-reppow   */
266   real       fudgeQQ;    /* The scaling factor for Coulomb 1-4: f*q1*q2  */
267   gmx_cmap_t cmap_grid;  /* The dihedral correction maps                 */
268 } gmx_ffparams_t;
269
270 enum {
271   ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
272 };
273
274 typedef struct
275 {
276   int ntypes;
277   int atnr;
278   t_functype *functype;
279   t_iparams  *iparams;
280   real fudgeQQ;
281   gmx_cmap_t cmap_grid;
282   t_iparams  *iparams_posres;
283   int iparams_posres_nalloc;
284
285   t_ilist il[F_NRE];
286   int ilsort;
287 } t_idef;
288
289 /*
290  * The struct t_idef defines all the interactions for the complete
291  * simulation. The structure is setup in such a way that the multinode
292  * version of the program  can use it as easy as the single node version.
293  * General field description:
294  *   int ntypes
295  *      defines the number of elements in functype[] and param[].
296  *   int nodeid
297  *      the node id (if parallel machines)
298  *   int atnr
299  *      the number of atomtypes
300  *   t_functype *functype
301  *      array of length ntypes, defines for every force type what type of 
302  *      function to use. Every "bond" with the same function but different 
303  *      force parameters is a different force type. The type identifier in the 
304  *      forceatoms[] array is an index in this array.
305  *   t_iparams *iparams
306  *      array of length ntypes, defines the parameters for every interaction
307  *      type. The type identifier in the actual interaction list
308  *      (ilist[ftype].iatoms[]) is an index in this array.
309  *   gmx_cmap_t cmap_grid
310  *      the grid for the dihedral pair correction maps.
311  *   t_iparams *iparams_posres
312  *      defines the parameters for position restraints only.
313  *      Position restraints are the only interactions that have different
314  *      parameters (reference positions) for different molecules
315  *      of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
316  *   t_ilist il[F_NRE]
317  *      The list of interactions for each type. Note that some,
318  *      such as LJ and COUL will have 0 entries.
319  */
320
321 typedef struct {
322   int  n;         /* n+1 is the number of points */
323   real scale;     /* distance between two points */
324   real *data;     /* the actual table data, per point there are 4 numbers */
325 } bondedtable_t;
326
327 #ifdef __cplusplus
328 }
329 #endif
330
331
332 #endif