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