Merge release-4-6 into master
[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 #include "nb_verlet.h"
41 #include "interaction_const.h"
42 #include "hw_info.h"
43
44 #ifdef __cplusplus
45 extern "C" {
46 #endif
47 #if 0
48 } /* fixes auto-indentation problems */
49 #endif
50
51 /* Abstract type for PME that is defined only in the routine that use them. */
52 typedef struct gmx_pme *gmx_pme_t;
53
54
55
56 /* Structure describing the data in a single table */
57 typedef struct
58 {
59     enum gmx_table_interaction  interaction; /* Types of interactions stored in this table */
60     enum gmx_table_format       format;      /* Interpolation type and data format */
61
62     real                        r;         /* range of the table */
63     int                         n;         /* n+1 is the number of table points */
64     real                        scale;     /* distance (nm) between two table points */
65     real                        scale_exp; /* distance for exponential part of VdW table, not always used */
66     real *                      data;      /* the actual table data */
67
68     /* Some information about the table layout. This can also be derived from the interpolation
69      * type and the table interactions, but it is convenient to have here for sanity checks, and it makes it
70      * much easier to access the tables in the nonbonded kernels when we can set the data from variables.
71      * It is always true that stride = formatsize*ninteractions
72      */
73     int                         formatsize;    /* Number of fp variables for each table point (1 for F, 2 for VF, 4 for YFGH, etc.) */
74     int                         ninteractions; /* Number of interactions in table, 1 for coul-only, 3 for coul+rep+disp. */
75     int                         stride;        /* Distance to next table point (number of fp variables per table point in total) */
76 } t_forcetable;
77
78 typedef struct
79 {
80     t_forcetable   table_elec;
81     t_forcetable   table_vdw;
82     t_forcetable   table_elec_vdw;
83
84     /* The actual neighbor lists, short and long range, see enum above
85      * for definition of neighborlist indices.
86      */
87     t_nblist nlist_sr[eNL_NR];
88     t_nblist nlist_lr[eNL_NR];
89 } t_nblists;
90
91 /* macros for the cginfo data in forcerec */
92 /* The maximum cg size in cginfo is 63
93  * because we only have space for 6 bits in cginfo,
94  * this cg size entry is actually only read with domain decomposition.
95  * But there is a smaller limit due to the t_excl data structure
96  * which is defined in nblist.h.
97  */
98 #define SET_CGINFO_GID(cgi,gid)      (cgi) = (((cgi)  &  ~65535)  |  (gid)   )
99 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   65535)
100 #define SET_CGINFO_EXCL_INTRA(cgi)   (cgi) =  ((cgi)  |  (1<<16))
101 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi)            &  (1<<16))
102 #define SET_CGINFO_EXCL_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<17))
103 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi)            &  (1<<17))
104 #define SET_CGINFO_SOLOPT(cgi,opt)   (cgi) = (((cgi)  & ~(3<<18)) | ((opt)<<18))
105 #define GET_CGINFO_SOLOPT(cgi)     (((cgi)>>18)       &   3)
106 #define SET_CGINFO_CONSTR(cgi)       (cgi) =  ((cgi)  |  (1<<20))
107 #define GET_CGINFO_CONSTR(cgi)     ( (cgi)            &  (1<<20))
108 #define SET_CGINFO_SETTLE(cgi)       (cgi) =  ((cgi)  |  (1<<21))
109 #define GET_CGINFO_SETTLE(cgi)     ( (cgi)            &  (1<<21))
110 /* This bit is only used with bBondComm in the domain decomposition */
111 #define SET_CGINFO_BOND_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<22))
112 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi)            &  (1<<22))
113 #define SET_CGINFO_HAS_VDW(cgi)      (cgi) =  ((cgi)  |  (1<<23))
114 #define GET_CGINFO_HAS_VDW(cgi)    ( (cgi)            &  (1<<23))
115 #define SET_CGINFO_HAS_Q(cgi)        (cgi) =  ((cgi)  |  (1<<24))
116 #define GET_CGINFO_HAS_Q(cgi)      ( (cgi)            &  (1<<24))
117 #define SET_CGINFO_NATOMS(cgi,opt)   (cgi) = (((cgi)  & ~(63<<25)) | ((opt)<<25))
118 #define GET_CGINFO_NATOMS(cgi)     (((cgi)>>25)       &   63)
119
120
121 /* Value to be used in mdrun for an infinite cut-off.
122  * Since we need to compare with the cut-off squared,
123  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
124  */
125 #define GMX_CUTOFF_INF 1E+18
126
127 /* enums for the neighborlist type */
128 enum { enbvdwNONE,enbvdwLJ,enbvdwBHAM,enbvdwTAB,enbvdwNR};
129 /* OOR is "one over r" -- standard coul */
130 enum { enbcoulNONE,enbcoulOOR,enbcoulRF,enbcoulTAB,enbcoulGB,enbcoulFEWALD,enbcoulNR};
131
132 enum { egCOULSR, egLJSR, egBHAMSR, egCOULLR, egLJLR, egBHAMLR,
133        egCOUL14, egLJ14, egGB, egNR };
134
135 typedef struct {
136   int  nener;        /* The number of energy group pairs     */
137   real *ener[egNR];  /* Energy terms for each pair of groups */
138 } gmx_grppairener_t;
139
140 typedef struct {
141   real term[F_NRE];    /* The energies for all different interaction types */
142   gmx_grppairener_t grpp;
143   double dvdl_lin[efptNR];       /* Contributions to dvdl with linear lam-dependence */
144   double dvdl_nonlin[efptNR];    /* Idem, but non-linear dependence                  */
145   int    n_lambda;
146   int    fep_state;              /*current fep state -- just for printing */
147   double *enerpart_lambda; /* Partial energy for lambda and flambda[] */
148 } gmx_enerdata_t;
149 /* The idea is that dvdl terms with linear lambda dependence will be added
150  * automatically to enerpart_lambda. Terms with non-linear lambda dependence
151  * should explicitly determine the energies at foreign lambda points
152  * when n_lambda > 0.
153  */
154
155 typedef struct {
156   int cg_start;
157   int cg_end;
158   int cg_mod;
159   int *cginfo;
160 } cginfo_mb_t;
161
162
163 /* ewald table type */
164 typedef struct ewald_tab *ewald_tab_t; 
165
166 typedef struct {
167     rvec *f;
168     int  f_nalloc;
169     unsigned red_mask; /* Mask for marking which parts of f are filled */
170     rvec *fshift;
171     real ener[F_NRE];
172     gmx_grppairener_t grpp;
173     real Vcorr;
174     real dvdl[efptNR];
175     tensor vir;
176 } f_thread_t;
177
178 typedef struct {
179   interaction_const_t *ic;
180
181   /* Domain Decomposition */
182   gmx_bool bDomDec;
183
184   /* PBC stuff */
185   int  ePBC;
186   gmx_bool bMolPBC;
187   int  rc_scaling;
188   rvec posres_com;
189   rvec posres_comB;
190
191   gmx_hw_info_t *hwinfo;
192   gmx_bool      use_cpu_acceleration;
193
194   /* Interaction for calculated in kernels. In many cases this is similar to
195    * the electrostatics settings in the inputrecord, but the difference is that
196    * these variables always specify the actual interaction in the kernel - if
197    * we are tabulating reaction-field the inputrec will say reaction-field, but
198    * the kernel interaction will say cubic-spline-table. To be safe we also
199    * have a kernel-specific setting for the modifiers - if the interaction is
200    * tabulated we already included the inputrec modification there, so the kernel
201    * modification setting will say 'none' in that case.
202    */
203   int nbkernel_elec_interaction;
204   int nbkernel_vdw_interaction;
205   int nbkernel_elec_modifier;
206   int nbkernel_vdw_modifier;
207
208   /* Use special N*N kernels? */
209   gmx_bool bAllvsAll;
210   /* Private work data */
211   void *AllvsAll_work;
212   void *AllvsAll_workgb;
213
214   /* Cut-Off stuff.
215    * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
216    */
217   real rlist,rlistlong;
218
219   /* Dielectric constant resp. multiplication factor for charges */
220   real zsquare,temp;
221   real epsilon_r,epsilon_rf,epsfac;  
222   
223   /* Constants for reaction fields */
224   real kappa,k_rf,c_rf;
225
226   /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
227   double qsum[2];
228   double q2sum[2];
229   rvec   mu_tot[2];
230
231   /* Dispersion correction stuff */
232   int  eDispCorr;
233
234   /* The shift of the shift or user potentials */
235   real enershiftsix;
236   real enershifttwelve;
237   /* Integrated differces for energy and virial with cut-off functions */
238   real enerdiffsix;
239   real enerdifftwelve;
240   real virdiffsix;
241   real virdifftwelve;
242   /* Constant for long range dispersion correction (average dispersion)
243    * for topology A/B ([0]/[1]) */
244   real avcsix[2];
245   /* Constant for long range repulsion term. Relative difference of about 
246    * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
247    */
248   real avctwelve[2];
249   
250   /* Fudge factors */
251   real fudgeQQ;
252
253   /* Table stuff */
254   gmx_bool bcoultab;
255   gmx_bool bvdwtab;
256   /* The normal tables are in the nblists struct(s) below */
257   t_forcetable tab14; /* for 1-4 interactions only */
258
259   /* PPPM & Shifting stuff */
260   int coulomb_modifier;
261   real rcoulomb_switch,rcoulomb;
262   real *phi;
263
264   /* VdW stuff */
265   int vdw_modifier;
266   double reppow;
267   real rvdw_switch,rvdw;
268   real bham_b_max;
269
270   /* Free energy */
271   int  efep;
272   real sc_alphavdw;
273   real sc_alphacoul;
274   int  sc_power;
275   real sc_r_power;
276   real sc_sigma6_def;
277   real sc_sigma6_min;
278   gmx_bool bSepDVDL;
279
280   /* NS Stuff */
281   int  eeltype;
282   int  vdwtype;
283   int  cg0,hcg;
284   /* solvent_opt contains the enum for the most common solvent
285    * in the system, which will be optimized.
286    * It can be set to esolNO to disable all water optimization */
287   int  solvent_opt;
288   int  nWatMol;
289   gmx_bool bGrid;
290   gmx_bool bExcl_IntraCGAll_InterCGNone;
291   cginfo_mb_t *cginfo_mb;
292   int  *cginfo;
293   rvec *cg_cm;
294   int  cg_nalloc;
295   rvec *shift_vec;
296
297   /* The neighborlists including tables */
298   int  nnblists;
299   int  *gid2nblists;
300   t_nblists *nblists;
301
302   int cutoff_scheme; /* group- or Verlet-style cutoff */
303   gmx_bool bNonbonded;    /* true if nonbonded calculations are *not* turned off */
304   nonbonded_verlet_t *nbv;
305
306   /* The wall tables (if used) */
307   int  nwall;
308   t_forcetable **wall_tab;
309
310   /* The number of charge groups participating in do_force_lowlevel */
311   int ncg_force;
312   /* The number of atoms participating in do_force_lowlevel */
313   int natoms_force;
314   /* The number of atoms participating in force and constraints */
315   int natoms_force_constr;
316   /* The allocation size of vectors of size natoms_force */
317   int nalloc_force;
318
319   /* Twin Range stuff, f_twin has size natoms_force */
320   gmx_bool bTwinRange;
321   int  nlr;
322   rvec *f_twin;
323
324   /* Forces that should not enter into the virial summation:
325    * PPPM/PME/Ewald/posres
326    */
327   gmx_bool bF_NoVirSum;
328   int  f_novirsum_n;
329   int  f_novirsum_nalloc;
330   rvec *f_novirsum_alloc;
331   /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
332    * points to the normal force vectors wen pressure is not requested.
333    */
334   rvec *f_novirsum;
335
336   /* Long-range forces and virial for PPPM/PME/Ewald */
337   gmx_pme_t pmedata;
338   tensor    vir_el_recip;
339
340   /* PME/Ewald stuff */
341   gmx_bool bEwald;
342   real ewaldcoeff;
343   ewald_tab_t ewald_table;
344
345   /* Virial Stuff */
346   rvec *fshift;
347   rvec vir_diag_posres;
348   dvec vir_wall_z;
349
350   /* Non bonded Parameter lists */
351   int  ntype; /* Number of atom types */
352   gmx_bool bBHAM;
353   real *nbfp;
354
355   /* Energy group pair flags */
356   int *egp_flags;
357
358   /* xmdrun flexible constraints */
359   real fc_stepsize;
360
361   /* Generalized born implicit solvent */
362   gmx_bool bGB;
363   /* Generalized born stuff */
364   real gb_epsilon_solvent;
365   /* Table data for GB */
366   t_forcetable gbtab;
367   /* VdW radius for each atomtype (dim is thus ntype) */
368   real *atype_radius;
369   /* Effective radius (derived from effective volume) for each type */
370   real *atype_vol;
371   /* Implicit solvent - surface tension for each atomtype */
372   real *atype_surftens;
373   /* Implicit solvent - radius for GB calculation */
374   real *atype_gb_radius;
375   /* Implicit solvent - overlap for HCT model */
376   real *atype_S_hct;
377   /* Generalized born interaction data */
378   gmx_genborn_t *born;
379
380   /* Table scale for GB */
381   real gbtabscale;
382   /* Table range for GB */
383   real gbtabr;
384   /* GB neighborlists (the sr list will contain for each atom all other atoms                                            
385    * (for use in the SA calculation) and the lr list will contain                                                  
386    * for each atom all atoms 1-4 or greater (for use in the GB calculation)                                        
387    */
388   t_nblist gblist_sr;
389   t_nblist gblist_lr;
390   t_nblist gblist;
391
392   /* Inverse square root of the Born radii for implicit solvent */
393   real *invsqrta;
394   /* Derivatives of the potential with respect to the Born radii */
395   real *dvda;
396   /* Derivatives of the Born radii with respect to coordinates */
397   real *dadx;
398   real *dadx_rawptr;
399   int   nalloc_dadx; /* Allocated size of dadx */
400         
401   /* If > 0 signals Test Particle Insertion,
402    * the value is the number of atoms of the molecule to insert
403    * Only the energy difference due to the addition of the last molecule
404    * should be calculated.
405    */
406   gmx_bool n_tpi;
407
408   /* Neighbor searching stuff */
409   gmx_ns_t ns;
410
411   /* QMMM stuff */
412   gmx_bool         bQMMM;
413   t_QMMMrec    *qr;
414
415   /* QM-MM neighborlists */
416   t_nblist QMMMlist;
417
418   /* Limit for printing large forces, negative is don't print */
419   real print_force;
420
421   /* coarse load balancing time measurement */
422   double t_fnbf;
423   double t_wait;
424   int timesteps;
425
426   /* parameter needed for AdResS simulation */
427   int  adress_type;
428   gmx_bool badress_tf_full_box;
429   real adress_const_wf;
430   real adress_ex_width;
431   real adress_hy_width;
432   int  adress_icor;
433   int  adress_site;
434   rvec adress_refs;
435   int n_adress_tf_grps;
436   int * adress_tf_table_index;
437   int *adress_group_explicit;
438   t_forcetable *  atf_tabs;
439   real adress_ex_forcecap;
440   gmx_bool adress_do_hybridpairs;
441
442   /* User determined parameters, copied from the inputrec */
443   int  userint1;
444   int  userint2;
445   int  userint3;
446   int  userint4;
447   real userreal1;
448   real userreal2;
449   real userreal3;
450   real userreal4;
451
452   /* Thread local force and energy data */ 
453   /* FIXME move to bonded_thread_data_t */
454   int  nthreads;
455   int  red_ashift;
456   int  red_nblock;
457   f_thread_t *f_t;
458
459   /* Exclusion load distribution over the threads */
460   int  *excl_load;
461 } t_forcerec;
462
463 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
464  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
465  * in the code, but beware if you are using these macros externally.
466  */
467 #define C6(nbfp,ntp,ai,aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
468 #define C12(nbfp,ntp,ai,aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
469 #define BHAMC(nbfp,ntp,ai,aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
470 #define BHAMA(nbfp,ntp,ai,aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
471 #define BHAMB(nbfp,ntp,ai,aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
472
473 #ifdef __cplusplus
474 }
475 #endif
476