Merge "Merge branch 'release-4-6'"
[alexxy/gromacs.git] / src / gromacs / legacyheaders / types / forcerec.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 #include "ns.h"
37 #include "genborn.h"
38 #include "qmmmrec.h"
39 #include "idef.h"
40
41 #ifdef __cplusplus
42 extern "C" {
43 #endif
44
45 /* Abstract type for PME that is defined only in the routine that use them. */
46 typedef struct gmx_pme *gmx_pme_t;
47
48 typedef struct {
49   real r;         /* range of the table */
50   int  n;         /* n+1 is the number of points */
51   real scale;     /* distance between two points */
52   real scale_exp; /* distance for exponential Buckingham table */
53   real *tab;      /* the actual tables, per point there are  4 numbers for
54                    * Coulomb, dispersion and repulsion (in total 12 numbers)
55                    */
56 } t_forcetable;
57
58 typedef struct {
59   t_forcetable tab;
60   /* We duplicate tables for cache optimization purposes */
61   real *coultab;      /* Coul only */
62   real *vdwtab;       /* Vdw only   */
63   /* The actual neighbor lists, short and long range, see enum above
64    * for definition of neighborlist indices.
65    */
66   t_nblist nlist_sr[eNL_NR];
67   t_nblist nlist_lr[eNL_NR];
68 } t_nblists;
69
70 /* macros for the cginfo data in forcerec */
71 /* The maximum cg size in cginfo is 255,
72  * because we only have space for 8 bits in cginfo,
73  * this cg size entry is actually only read with domain decomposition.
74  * But there is a smaller limit due to the t_excl data structure
75  * which is defined in nblist.h.
76  */
77 #define SET_CGINFO_GID(cgi,gid)      (cgi) = (((cgi)  &  ~65535)  |  (gid)   )
78 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   65535)
79 #define SET_CGINFO_EXCL_INTRA(cgi)   (cgi) =  ((cgi)  |  (1<<16))
80 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi)            &  (1<<16))
81 #define SET_CGINFO_EXCL_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<17))
82 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi)            &  (1<<17))
83 #define SET_CGINFO_SOLOPT(cgi,opt)   (cgi) = (((cgi)  & ~(15<<18)) | ((opt)<<18))
84 #define GET_CGINFO_SOLOPT(cgi)     (((cgi)>>18)       &   15)
85 /* This bit is only used with bBondComm in the domain decomposition */
86 #define SET_CGINFO_BOND_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<22))
87 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi)            &  (1<<22))
88 #define SET_CGINFO_NATOMS(cgi,opt)   (cgi) = (((cgi)  & ~(255<<23)) | ((opt)<<23))
89 #define GET_CGINFO_NATOMS(cgi)     (((cgi)>>23)       &   255)
90
91
92 /* Value to be used in mdrun for an infinite cut-off.
93  * Since we need to compare with the cut-off squared,
94  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
95  */
96 #define GMX_CUTOFF_INF 1E+18
97
98 /* enums for the neighborlist type */
99 enum { enbvdwNONE,enbvdwLJ,enbvdwBHAM,enbvdwTAB,enbvdwNR};
100 /* OOR is "one over r" -- standard coul */
101 enum { enbcoulNONE,enbcoulOOR,enbcoulRF,enbcoulTAB,enbcoulGB,enbcoulFEWALD,enbcoulNR};
102
103 enum { egCOULSR, egLJSR, egBHAMSR, egCOULLR, egLJLR, egBHAMLR,
104        egCOUL14, egLJ14, egGB, egNR };
105
106 typedef struct {
107   int  nener;        /* The number of energy group pairs     */
108   real *ener[egNR];  /* Energy terms for each pair of groups */
109 } gmx_grppairener_t;
110
111 typedef struct {
112   real term[F_NRE];    /* The energies for all different interaction types */
113   gmx_grppairener_t grpp;
114   double dvdl_lin[efptNR];       /* Contributions to dvdl with linear lam-dependence */
115   double dvdl_nonlin[efptNR];    /* Idem, but non-linear dependence                  */
116   int    n_lambda;
117   int    fep_state;              /*current fep state -- just for printing */
118   double *enerpart_lambda; /* Partial energy for lambda and flambda[] */
119 } gmx_enerdata_t;
120 /* The idea is that dvdl terms with linear lambda dependence will be added
121  * automatically to enerpart_lambda. Terms with non-linear lambda dependence
122  * should explicitly determine the energies at foreign lambda points
123  * when n_lambda > 0.
124  */
125
126 typedef struct {
127   int cg_start;
128   int cg_end;
129   int cg_mod;
130   int *cginfo;
131 } cginfo_mb_t;
132
133
134 /* ewald table type */
135 typedef struct ewald_tab *ewald_tab_t; 
136
137 typedef struct {
138   /* Domain Decomposition */
139   gmx_bool bDomDec;
140
141   /* PBC stuff */
142   int  ePBC;
143   gmx_bool bMolPBC;
144   int  rc_scaling;
145   rvec posres_com;
146   rvec posres_comB;
147
148   gmx_bool UseOptimizedKernels;
149
150   /* Use special N*N kernels? */
151   gmx_bool bAllvsAll;
152   /* Private work data */
153   void *AllvsAll_work;
154   void *AllvsAll_workgb;
155
156   /* Cut-Off stuff.
157    * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
158    */
159   real rlist,rlistlong;
160   
161   /* Dielectric constant resp. multiplication factor for charges */
162   real zsquare,temp;
163   real epsilon_r,epsilon_rf,epsfac;  
164   
165   /* Constants for reaction fields */
166   real kappa,k_rf,c_rf;
167
168   /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
169   double qsum[2];
170   rvec   mu_tot[2];
171
172   /* Dispersion correction stuff */
173   int  eDispCorr;
174   /* The shift of the shift or user potentials */
175   real enershiftsix;
176   real enershifttwelve;
177   /* Integrated differces for energy and virial with cut-off functions */
178   real enerdiffsix;
179   real enerdifftwelve;
180   real virdiffsix;
181   real virdifftwelve;
182   /* Constant for long range dispersion correction (average dispersion)
183    * for topology A/B ([0]/[1]) */
184   real avcsix[2];
185   /* Constant for long range repulsion term. Relative difference of about 
186    * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
187    */
188   real avctwelve[2];
189   
190   /* Fudge factors */
191   real fudgeQQ;
192
193   /* Table stuff */
194   gmx_bool bcoultab;
195   gmx_bool bvdwtab;
196   /* The normal tables are in the nblists struct(s) below */
197   t_forcetable tab14; /* for 1-4 interactions only */
198
199   /* PPPM & Shifting stuff */
200   real rcoulomb_switch,rcoulomb;
201   real *phi;
202
203   /* VdW stuff */
204   double reppow;
205   real rvdw_switch,rvdw;
206   real bham_b_max;
207
208   /* Free energy */
209   int  efep;
210   real sc_alphavdw;
211   real sc_alphacoul;
212   int  sc_power;
213   real sc_r_power;
214   real sc_sigma6_def;
215   real sc_sigma6_min;
216   gmx_bool bSepDVDL;
217
218   /* NS Stuff */
219   int  eeltype;
220   int  vdwtype;
221   int  cg0,hcg;
222   /* solvent_opt contains the enum for the most common solvent
223    * in the system, which will be optimized.
224    * It can be set to esolNO to disable all water optimization */
225   int  solvent_opt;
226   int  nWatMol;
227   gmx_bool bGrid;
228   cginfo_mb_t *cginfo_mb;
229   int  *cginfo;
230   rvec *cg_cm;
231   int  cg_nalloc;
232   rvec *shift_vec;
233
234   /* The neighborlists including tables */
235   int  nnblists;
236   int  *gid2nblists;
237   t_nblists *nblists;
238
239   /* The wall tables (if used) */
240   int  nwall;
241   t_forcetable **wall_tab;
242
243   /* This mask array of length nn determines whether or not this bit of the
244    * neighbourlists should be computed. Usually all these are true of course,
245    * but not when shells are used. During minimisation all the forces that 
246    * include shells are done, then after minimsation is converged the remaining
247    * forces are computed.
248    */
249   /* gmx_bool *bMask; */
250
251   /* The number of charge groups participating in do_force_lowlevel */
252   int ncg_force;
253   /* The number of atoms participating in do_force_lowlevel */
254   int natoms_force;
255   /* The number of atoms participating in force and constraints */
256   int natoms_force_constr;
257   /* The allocation size of vectors of size natoms_force */
258   int nalloc_force;
259
260   /* Twin Range stuff, f_twin has size natoms_force */
261   gmx_bool bTwinRange;
262   int  nlr;
263   rvec *f_twin;
264
265   /* Forces that should not enter into the virial summation:
266    * PPPM/PME/Ewald/posres
267    */
268   gmx_bool bF_NoVirSum;
269   int  f_novirsum_n;
270   int  f_novirsum_nalloc;
271   rvec *f_novirsum_alloc;
272   /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
273    * points to the normal force vectors wen pressure is not requested.
274    */
275   rvec *f_novirsum;
276
277   /* Long-range forces and virial for PPPM/PME/Ewald */
278   gmx_pme_t pmedata;
279   tensor    vir_el_recip;
280
281   /* PME/Ewald stuff */
282   gmx_bool bEwald;
283   real ewaldcoeff;
284   ewald_tab_t ewald_table;
285
286   /* Virial Stuff */
287   rvec *fshift;
288   rvec vir_diag_posres;
289   dvec vir_wall_z;
290
291   /* Non bonded Parameter lists */
292   int  ntype; /* Number of atom types */
293   gmx_bool bBHAM;
294   real *nbfp;
295
296   /* Energy group pair flags */
297   int *egp_flags;
298
299   /* xmdrun flexible constraints */
300   real fc_stepsize;
301
302   /* Generalized born implicit solvent */
303   gmx_bool bGB;
304   /* Generalized born stuff */
305   real gb_epsilon_solvent;
306   /* Table data for GB */
307   t_forcetable gbtab;
308   /* VdW radius for each atomtype (dim is thus ntype) */
309   real *atype_radius;
310   /* Effective radius (derived from effective volume) for each type */
311   real *atype_vol;
312   /* Implicit solvent - surface tension for each atomtype */
313   real *atype_surftens;
314   /* Implicit solvent - radius for GB calculation */
315   real *atype_gb_radius;
316   /* Implicit solvent - overlap for HCT model */
317   real *atype_S_hct;
318   /* Generalized born interaction data */
319   gmx_genborn_t *born;
320
321   /* Table scale for GB */
322   real gbtabscale;
323   /* Table range for GB */
324   real gbtabr;
325   /* GB neighborlists (the sr list will contain for each atom all other atoms                                            
326    * (for use in the SA calculation) and the lr list will contain                                                  
327    * for each atom all atoms 1-4 or greater (for use in the GB calculation)                                        
328    */
329   t_nblist gblist_sr;
330   t_nblist gblist_lr;
331   t_nblist gblist;
332
333   /* Inverse square root of the Born radii for implicit solvent */
334   real *invsqrta;
335   /* Derivatives of the potential with respect to the Born radii */
336   real *dvda;
337   /* Derivatives of the Born radii with respect to coordinates */
338   real *dadx;
339   real *dadx_rawptr;
340   int   nalloc_dadx; /* Allocated size of dadx */
341         
342   /* If > 0 signals Test Particle Insertion,
343    * the value is the number of atoms of the molecule to insert
344    * Only the energy difference due to the addition of the last molecule
345    * should be calculated.
346    */
347   gmx_bool n_tpi;
348
349   /* Neighbor searching stuff */
350   gmx_ns_t ns;
351
352   /* QMMM stuff */
353   gmx_bool         bQMMM;
354   t_QMMMrec    *qr;
355
356   /* QM-MM neighborlists */
357   t_nblist QMMMlist;
358
359   /* Limit for printing large forces, negative is don't print */
360   real print_force;
361
362   /* coarse load balancing time measurement */
363   double t_fnbf;
364   double t_wait;
365   int timesteps;
366
367   /* parameter needed for AdResS simulation */
368   int  adress_type;
369   gmx_bool badress_tf_full_box;
370   real adress_const_wf;
371   real adress_ex_width;
372   real adress_hy_width;
373   int  adress_icor;
374   int  adress_site;
375   rvec adress_refs;
376   int n_adress_tf_grps;
377   int * adress_tf_table_index;
378   int *adress_group_explicit;
379   t_forcetable *  atf_tabs;
380   real adress_ex_forcecap;
381   gmx_bool adress_do_hybridpairs;
382
383   /* User determined parameters, copied from the inputrec */
384   int  userint1;
385   int  userint2;
386   int  userint3;
387   int  userint4;
388   real userreal1;
389   real userreal2;
390   real userreal3;
391   real userreal4;
392 } t_forcerec;
393
394 #define C6(nbfp,ntp,ai,aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
395 #define C12(nbfp,ntp,ai,aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
396 #define BHAMC(nbfp,ntp,ai,aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
397 #define BHAMA(nbfp,ntp,ai,aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
398 #define BHAMB(nbfp,ntp,ai,aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
399
400 #ifdef __cplusplus
401 }
402 #endif
403