Added two new bonded functions.
[alexxy/gromacs.git] / include / 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_DISRES,
108   F_DISRESVIOL,
109   F_ORIRES,
110   F_ORIRESDEV,
111   F_ANGRES,
112   F_ANGRESZ,
113   F_DIHRES,
114   F_DIHRESVIOL,
115   F_CONSTR,
116   F_CONSTRNC,
117   F_SETTLE,
118   F_VSITE2,
119   F_VSITE3,
120   F_VSITE3FD,
121   F_VSITE3FAD,
122   F_VSITE3OUT,
123   F_VSITE4FD,
124   F_VSITE4FDN,
125   F_VSITEN,
126   F_COM_PULL,
127   F_EQM,
128   F_EPOT,
129   F_EKIN,
130   F_ETOT,
131   F_ECONSERVED,
132   F_TEMP,
133   F_VTEMP,
134   F_PDISPCORR,
135   F_PRES,
136   F_DVDL,
137   F_DKDL,
138   F_DHDL_CON,
139   F_NRE         /* This number is for the total number of energies      */
140 };
141   
142 typedef union
143 {
144   /* Some parameters have A and B values for free energy calculations.
145    * The B values are not used for regular simulations of course.
146    * Free Energy for nonbondeds can be computed by changing the atom type.
147    * The harmonic type is used for all harmonic potentials:
148    * bonds, angles and improper dihedrals
149    */
150   struct {real a,b,c;                                      } bham;
151   struct {real rA,krA,rB,krB;                              } harmonic;
152   struct {real klinA,aA,klinB,aB;                          } linangle;
153   struct {real lowA,up1A,up2A,kA,lowB,up1B,up2B,kB;        } restraint;
154   /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */ 
155   struct {real b0,kb,kcub;                                 } cubic;
156   struct {real bm,kb;                                      } fene;
157   struct {real r1e,r2e,krr;                                } cross_bb;
158   struct {real r1e,r2e,r3e,krt;                            } cross_ba;
159   struct {real theta,ktheta,r13,kUB;                       } u_b;
160   struct {real theta,c[5];                                 } qangle; 
161   struct {real alpha;                                      } polarize;
162   struct {real alpha,drcut,khyp;                           } anharm_polarize;
163   struct {real al_x,al_y,al_z,rOH,rHH,rOD;                 } wpol;
164   struct {real a,alpha1,alpha2,rfac;                       } thole;
165   struct {real c6,c12;                                     } lj;
166   struct {real c6A,c12A,c6B,c12B;                          } lj14;
167   struct {real fqq,qi,qj,c6,c12;                           } ljc14;
168   struct {real qi,qj,c6,c12;                               } ljcnb;
169   /* Proper dihedrals can not have different multiplicity when
170    * doing free energy calculations, because the potential would not
171    * be periodic anymore.
172    */ 
173   struct {real phiA,cpA;int mult;real phiB,cpB;            } pdihs;
174   struct {real dA,dB;                                      } constr;
175   /* Settle can not be used for Free energy calculations of water bond geometry.
176    * Use shake (or lincs) instead if you have to change the water bonds.
177    */
178   struct {real doh,dhh;                                   } settle;
179   /* No free energy supported for morse bonds */ 
180   struct {real b0,cb,beta;                                } morse;
181   struct {real pos0A[DIM],fcA[DIM],pos0B[DIM],fcB[DIM];   } posres;
182   struct {real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];          } rbdihs;
183   struct {real a,b,c,d,e,f;                               } vsite;   
184   struct {int  n; real a;                                 } vsiten;   
185   /* NOTE: npair is only set after reading the tpx file */
186   struct {real low,up1,up2,kfac;int type,label,npair;     } disres; 
187   struct {real phi,dphi,kfac;int label,power;             } dihres;  
188   struct {int  ex,power,label; real c,obs,kfac;           } orires;
189   struct {int  table;real kA;real kB;                     } tab;
190   struct {real sar,st,pi,gbr,bmlt;                        } gb;
191   struct {int cmapA,cmapB;                                } cmap;
192   struct {real buf[MAXFORCEPARAM];                        } generic; /* Conversion */
193 } t_iparams;
194
195 typedef int t_functype;
196
197 /*
198  * The nonperturbed/perturbed interactions are now separated (sorted) in the
199  * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and 
200  * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded 
201  * interactions.
202  */
203 typedef struct
204 {
205   int nr;
206   int nr_nonperturbed;
207   t_iatom *iatoms;
208   int nalloc;
209 } t_ilist;
210
211 /*
212  * The struct t_ilist defines a list of atoms with their interactions. 
213  * General field description:
214  *   int nr
215  *      the size (nr elements) of the interactions array (iatoms[]).
216  *   t_iatom *iatoms
217  *      specifies which atoms are involved in an interaction of a certain 
218  *       type. The layout of this array is as follows:
219  *
220  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
221  *        |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
222  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
223  *
224  *      So for interaction type type1 3 atoms are needed, and for type2 and 
225  *      type3 only 2. The type identifier is used to select the function to 
226  *      calculate the interaction and its actual parameters. This type 
227  *      identifier is an index in a params[] and functype[] array.
228  */
229
230 typedef struct
231 {
232         real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
233         /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
234 } cmapdata_t;
235
236 typedef struct
237 {
238         int ngrid;            /* Number of allocated cmap (cmapdata_t ) grids */
239         int grid_spacing;     /* Grid spacing */
240         cmapdata_t *cmapdata; /* Pointer to grid with actual, pre-interpolated data */
241 } gmx_cmap_t;
242
243
244 typedef struct
245 {
246   int        ntypes;
247   int        atnr;
248   t_functype *functype;
249   t_iparams  *iparams;
250   double     reppow;     /* The repulsion power for VdW: C12*r^-reppow   */
251   real       fudgeQQ;    /* The scaling factor for Coulomb 1-4: f*q1*q2  */
252   gmx_cmap_t cmap_grid;  /* The dihedral correction maps                 */
253 } gmx_ffparams_t;
254
255 enum {
256   ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
257 };
258
259 typedef struct
260 {
261   int ntypes;
262   int atnr;
263   t_functype *functype;
264   t_iparams  *iparams;
265   real fudgeQQ;
266   gmx_cmap_t cmap_grid;
267   t_iparams  *iparams_posres;
268   int iparams_posres_nalloc;
269
270   t_ilist il[F_NRE];
271   int ilsort;
272 } t_idef;
273
274 /*
275  * The struct t_idef defines all the interactions for the complete
276  * simulation. The structure is setup in such a way that the multinode
277  * version of the program  can use it as easy as the single node version.
278  * General field description:
279  *   int ntypes
280  *      defines the number of elements in functype[] and param[].
281  *   int nodeid
282  *      the node id (if parallel machines)
283  *   int atnr
284  *      the number of atomtypes
285  *   t_functype *functype
286  *      array of length ntypes, defines for every force type what type of 
287  *      function to use. Every "bond" with the same function but different 
288  *      force parameters is a different force type. The type identifier in the 
289  *      forceatoms[] array is an index in this array.
290  *   t_iparams *iparams
291  *      array of length ntypes, defines the parameters for every interaction
292  *      type. The type identifier in the actual interaction list
293  *      (ilist[ftype].iatoms[]) is an index in this array.
294  *   gmx_cmap_t cmap_grid
295  *      the grid for the dihedral pair correction maps.
296  *   t_iparams *iparams_posres
297  *      defines the parameters for position restraints only.
298  *      Position restraints are the only interactions that have different
299  *      parameters (reference positions) for different molecules
300  *      of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
301  *   t_ilist il[F_NRE]
302  *      The list of interactions for each type. Note that some,
303  *      such as LJ and COUL will have 0 entries.
304  */
305
306 typedef struct {
307   int  n;         /* n+1 is the number of points */
308   real scale;     /* distance between two points */
309   real *tab;      /* the actual tables, per point there are  4 numbers */
310 } bondedtable_t;
311
312 #ifdef __cplusplus
313 }
314 #endif
315
316
317 #endif