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