Finish removing group-scheme all-vs-all support
[alexxy/gromacs.git] / src / gromacs / mdtypes / 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,2016,2017,2018,2019, 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 #ifndef GMX_MDTYPES_TYPES_FORCEREC_H
38 #define GMX_MDTYPES_TYPES_FORCEREC_H
39
40 #include <array>
41 #include <vector>
42
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdtypes/interaction_const.h"
45 #include "gromacs/mdtypes/md_enums.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/real.h"
48
49 struct ForceProviders;
50
51 /* Abstract type for PME that is defined only in the routine that use them. */
52 struct gmx_ns_t;
53 struct gmx_pme_t;
54 struct nonbonded_verlet_t;
55 struct bonded_threading_t;
56 struct t_forcetable;
57 struct t_nblist;
58 struct t_nblists;
59 struct t_QMMMrec;
60
61 namespace gmx
62 {
63 class GpuBonded;
64 }
65
66 /* macros for the cginfo data in forcerec
67  *
68  * Since the tpx format support max 256 energy groups, we do the same here.
69  * Note that we thus have bits 8-14 still unused.
70  *
71  * The maximum cg size in cginfo is 63
72  * because we only have space for 6 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)  &  ~255) | (gid))
78 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   255)
79 #define SET_CGINFO_FEP(cgi)          (cgi) =  ((cgi)  |  (1<<15))
80 #define GET_CGINFO_FEP(cgi)        ( (cgi)            &  (1<<15))
81 #define SET_CGINFO_EXCL_INTRA(cgi)   (cgi) =  ((cgi)  |  (1<<16))
82 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi)            &  (1<<16))
83 #define SET_CGINFO_EXCL_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<17))
84 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi)            &  (1<<17))
85 #define SET_CGINFO_SOLOPT(cgi, opt)  (cgi) = (((cgi)  & ~(3<<18)) | ((opt)<<18))
86 #define GET_CGINFO_SOLOPT(cgi)     (((cgi)>>18)       &   3)
87 #define SET_CGINFO_CONSTR(cgi)       (cgi) =  ((cgi)  |  (1<<20))
88 #define GET_CGINFO_CONSTR(cgi)     ( (cgi)            &  (1<<20))
89 #define SET_CGINFO_SETTLE(cgi)       (cgi) =  ((cgi)  |  (1<<21))
90 #define GET_CGINFO_SETTLE(cgi)     ( (cgi)            &  (1<<21))
91 /* This bit is only used with bBondComm in the domain decomposition */
92 #define SET_CGINFO_BOND_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<22))
93 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi)            &  (1<<22))
94 #define SET_CGINFO_HAS_VDW(cgi)      (cgi) =  ((cgi)  |  (1<<23))
95 #define GET_CGINFO_HAS_VDW(cgi)    ( (cgi)            &  (1<<23))
96 #define SET_CGINFO_HAS_Q(cgi)        (cgi) =  ((cgi)  |  (1<<24))
97 #define GET_CGINFO_HAS_Q(cgi)      ( (cgi)            &  (1<<24))
98 #define SET_CGINFO_NATOMS(cgi, opt)  (cgi) = (((cgi)  & ~(63<<25)) | ((opt)<<25))
99 #define GET_CGINFO_NATOMS(cgi)     (((cgi)>>25)       &   63)
100
101
102 /* Value to be used in mdrun for an infinite cut-off.
103  * Since we need to compare with the cut-off squared,
104  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
105  */
106 #define GMX_CUTOFF_INF 1E+18
107
108 /* enums for the neighborlist type */
109 enum {
110     enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
111 };
112
113 struct cginfo_mb_t
114 {
115     int  cg_start;
116     int  cg_end;
117     int  cg_mod;
118     int *cginfo;
119 };
120
121
122 /* Forward declaration of type for managing Ewald tables */
123 struct gmx_ewald_tab_t;
124
125 struct ewald_corr_thread_t;
126
127 struct t_forcerec { // NOLINT (clang-analyzer-optin.performance.Padding)
128     struct interaction_const_t *ic;
129
130     /* PBC stuff */
131     int                         ePBC;
132     //! Whether PBC must be considered for bonded interactions.
133     gmx_bool                    bMolPBC;
134     int                         rc_scaling;
135     rvec                        posres_com;
136     rvec                        posres_comB;
137
138     gmx_bool                    use_simd_kernels;
139
140     /* Interaction for calculated in kernels. In many cases this is similar to
141      * the electrostatics settings in the inputrecord, but the difference is that
142      * these variables always specify the actual interaction in the kernel - if
143      * we are tabulating reaction-field the inputrec will say reaction-field, but
144      * the kernel interaction will say cubic-spline-table. To be safe we also
145      * have a kernel-specific setting for the modifiers - if the interaction is
146      * tabulated we already included the inputrec modification there, so the kernel
147      * modification setting will say 'none' in that case.
148      */
149     int nbkernel_elec_interaction;
150     int nbkernel_vdw_interaction;
151     int nbkernel_elec_modifier;
152     int nbkernel_vdw_modifier;
153
154     /* Cut-Off stuff.
155      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
156      */
157     real rlist;
158
159     /* Parameters for generalized reaction field */
160     real zsquare, temp;
161
162     /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
163     double qsum[2];
164     double q2sum[2];
165     double c6sum[2];
166     rvec   mu_tot[2];
167
168     /* Dispersion correction stuff */
169     int                  eDispCorr;
170     int                  numAtomsForDispersionCorrection;
171     struct t_forcetable *dispersionCorrectionTable;
172
173     /* The shift of the shift or user potentials */
174     real enershiftsix;
175     real enershifttwelve;
176     /* Integrated differces for energy and virial with cut-off functions */
177     real enerdiffsix;
178     real enerdifftwelve;
179     real virdiffsix;
180     real virdifftwelve;
181     /* Constant for long range dispersion correction (average dispersion)
182      * for topology A/B ([0]/[1]) */
183     real avcsix[2];
184     /* Constant for long range repulsion term. Relative difference of about
185      * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
186      */
187     real avctwelve[2];
188
189     /* Fudge factors */
190     real fudgeQQ;
191
192     /* Table stuff */
193     gmx_bool             bcoultab;
194     gmx_bool             bvdwtab;
195     /* The normal tables are in the nblists struct(s) below */
196
197     struct t_forcetable *pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
198
199     /* Free energy */
200     int      efep;
201     real     sc_alphavdw;
202     real     sc_alphacoul;
203     int      sc_power;
204     real     sc_r_power;
205     real     sc_sigma6_def;
206     real     sc_sigma6_min;
207
208     /* NS Stuff */
209     int  cg0, hcg;
210     /* solvent_opt contains the enum for the most common solvent
211      * in the system, which will be optimized.
212      * It can be set to esolNO to disable all water optimization */
213     int                 solvent_opt;
214     int                 nWatMol;
215     gmx_bool            bGrid;
216     gmx_bool            bExcl_IntraCGAll_InterCGNone;
217     struct cginfo_mb_t *cginfo_mb;
218     int                *cginfo;
219     rvec               *cg_cm;
220     int                 cg_nalloc;
221     rvec               *shift_vec;
222
223     /* The neighborlists including tables */
224     int                        nnblists;
225     int                       *gid2nblists;
226     struct t_nblists          *nblists;
227
228     int                        cutoff_scheme; /* group- or Verlet-style cutoff */
229     gmx_bool                   bNonbonded;    /* true if nonbonded calculations are *not* turned off */
230     struct nonbonded_verlet_t *nbv;
231
232     /* The wall tables (if used) */
233     int                    nwall;
234     struct t_forcetable ***wall_tab;
235
236     /* The number of charge groups participating in do_force_lowlevel */
237     int ncg_force;
238     /* The number of atoms participating in do_force_lowlevel */
239     int natoms_force;
240     /* The number of atoms participating in force and constraints */
241     int natoms_force_constr;
242     /* The allocation size of vectors of size natoms_force */
243     int nalloc_force;
244
245     /* Forces that should not enter into the coord x force virial summation:
246      * PPPM/PME/Ewald/posres/ForceProviders
247      */
248     /* True when we have contributions that are directly added to the virial */
249     gmx_bool                 haveDirectVirialContributions;
250     /* TODO: Replace the pointer by an object once we got rid of C */
251     std::vector<gmx::RVec>  *forceBufferForDirectVirialContributions;
252
253     /* Data for PPPM/PME/Ewald */
254     struct gmx_pme_t *pmedata;
255     int               ljpme_combination_rule;
256
257     /* PME/Ewald stuff */
258     struct gmx_ewald_tab_t *ewald_table;
259
260     /* Shift force array for computing the virial */
261     rvec *fshift;
262
263     /* Non bonded Parameter lists */
264     int      ntype; /* Number of atom types */
265     gmx_bool bBHAM;
266     real    *nbfp;
267     real    *ljpme_c6grid; /* C6-values used on grid in LJPME */
268
269     /* Energy group pair flags */
270     int *egp_flags;
271
272     /* Shell molecular dynamics flexible constraints */
273     real fc_stepsize;
274
275     /* If > 0 signals Test Particle Insertion,
276      * the value is the number of atoms of the molecule to insert
277      * Only the energy difference due to the addition of the last molecule
278      * should be calculated.
279      */
280     int n_tpi;
281
282     /* Neighbor searching stuff */
283     struct gmx_ns_t *ns;
284
285     /* QMMM stuff */
286     gmx_bool          bQMMM;
287     struct t_QMMMrec *qr;
288
289     /* QM-MM neighborlists */
290     struct t_nblist        *QMMMlist;
291
292     /* Limit for printing large forces, negative is don't print */
293     real print_force;
294
295     /* coarse load balancing time measurement */
296     double t_fnbf;
297     double t_wait;
298     int    timesteps;
299
300     /* User determined parameters, copied from the inputrec */
301     int  userint1;
302     int  userint2;
303     int  userint3;
304     int  userint4;
305     real userreal1;
306     real userreal2;
307     real userreal3;
308     real userreal4;
309
310     /* Pointer to struct for managing threading of bonded force calculation */
311     struct bonded_threading_t *bondedThreading;
312
313     /* TODO: Replace the pointer by an object once we got rid of C */
314     gmx::GpuBonded *gpuBonded;
315
316     /* Ewald correction thread local virial and energy data */
317     int                         nthread_ewc;
318     struct ewald_corr_thread_t *ewc_t;
319
320     struct ForceProviders      *forceProviders;
321 };
322
323 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
324  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
325  * in the code, but beware if you are using these macros externally.
326  */
327 #define C6(nbfp, ntp, ai, aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
328 #define C12(nbfp, ntp, ai, aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
329 #define BHAMC(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
330 #define BHAMA(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
331 #define BHAMB(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
332
333 #endif