ec01003d243b5308c65a4d29f91202b8d11075e3
[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 <memory>
42 #include <vector>
43
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/mdtypes/interaction_const.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
49
50 /* Abstract type for PME that is defined only in the routine that use them. */
51 struct gmx_ns_t;
52 struct gmx_pme_t;
53 struct nonbonded_verlet_t;
54 struct bonded_threading_t;
55 class DispersionCorrection;
56 struct t_forcetable;
57 struct t_QMMMrec;
58
59 namespace gmx
60 {
61 class GpuBonded;
62 class ForceProviders;
63 }
64
65 /* macros for the cginfo data in forcerec
66  *
67  * Since the tpx format support max 256 energy groups, we do the same here.
68  * Note that we thus have bits 8-14 still unused.
69  *
70  * The maximum cg size in cginfo is 63
71  * because we only have space for 6 bits in cginfo,
72  * this cg size entry is actually only read with domain decomposition.
73  */
74 #define SET_CGINFO_GID(cgi, gid)     (cgi) = (((cgi)  &  ~255) | (gid))
75 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   255)
76 #define SET_CGINFO_FEP(cgi)          (cgi) =  ((cgi)  |  (1<<15))
77 #define GET_CGINFO_FEP(cgi)        ( (cgi)            &  (1<<15))
78 #define SET_CGINFO_EXCL_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<17))
79 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi)            &  (1<<17))
80 #define SET_CGINFO_CONSTR(cgi)       (cgi) =  ((cgi)  |  (1<<20))
81 #define GET_CGINFO_CONSTR(cgi)     ( (cgi)            &  (1<<20))
82 #define SET_CGINFO_SETTLE(cgi)       (cgi) =  ((cgi)  |  (1<<21))
83 #define GET_CGINFO_SETTLE(cgi)     ( (cgi)            &  (1<<21))
84 /* This bit is only used with bBondComm in the domain decomposition */
85 #define SET_CGINFO_BOND_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<22))
86 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi)            &  (1<<22))
87 #define SET_CGINFO_HAS_VDW(cgi)      (cgi) =  ((cgi)  |  (1<<23))
88 #define GET_CGINFO_HAS_VDW(cgi)    ( (cgi)            &  (1<<23))
89 #define SET_CGINFO_HAS_Q(cgi)        (cgi) =  ((cgi)  |  (1<<24))
90 #define GET_CGINFO_HAS_Q(cgi)      ( (cgi)            &  (1<<24))
91
92
93 /* Value to be used in mdrun for an infinite cut-off.
94  * Since we need to compare with the cut-off squared,
95  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
96  */
97 #define GMX_CUTOFF_INF 1E+18
98
99 /* enums for the neighborlist type */
100 enum {
101     enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
102 };
103
104 struct cginfo_mb_t
105 {
106     int  cg_start;
107     int  cg_end;
108     int  cg_mod;
109     int *cginfo;
110 };
111
112
113 /* Forward declaration of type for managing Ewald tables */
114 struct gmx_ewald_tab_t;
115
116 struct ewald_corr_thread_t;
117
118 struct t_forcerec { // NOLINT (clang-analyzer-optin.performance.Padding)
119     struct interaction_const_t *ic = nullptr;
120
121     /* PBC stuff */
122     int                         ePBC = 0;
123     //! Tells whether atoms inside a molecule can be in different periodic images,
124     //  i.e. whether we need to take into account PBC when computing distances inside molecules.
125     //  This determines whether PBC must be considered for e.g. bonded interactions.
126     gmx_bool                    bMolPBC     = FALSE;
127     int                         rc_scaling  = 0;
128     rvec                        posres_com  = { 0 };
129     rvec                        posres_comB = { 0 };
130
131     gmx_bool                    use_simd_kernels = FALSE;
132
133     /* Interaction for calculated in kernels. In many cases this is similar to
134      * the electrostatics settings in the inputrecord, but the difference is that
135      * these variables always specify the actual interaction in the kernel - if
136      * we are tabulating reaction-field the inputrec will say reaction-field, but
137      * the kernel interaction will say cubic-spline-table. To be safe we also
138      * have a kernel-specific setting for the modifiers - if the interaction is
139      * tabulated we already included the inputrec modification there, so the kernel
140      * modification setting will say 'none' in that case.
141      */
142     int nbkernel_elec_interaction = 0;
143     int nbkernel_vdw_interaction  = 0;
144     int nbkernel_elec_modifier    = 0;
145     int nbkernel_vdw_modifier     = 0;
146
147     /* Cut-Off stuff.
148      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
149      */
150     real rlist = 0;
151
152     /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
153     double qsum[2]   = { 0 };
154     double q2sum[2]  = { 0 };
155     double c6sum[2]  = { 0 };
156     rvec   mu_tot[2] = { { 0 } };
157
158     /* Dispersion correction stuff */
159     std::unique_ptr<DispersionCorrection> dispersionCorrection;
160
161     /* Fudge factors */
162     real fudgeQQ = 0;
163
164     /* Table stuff */
165     gmx_bool             bcoultab = FALSE;
166     gmx_bool             bvdwtab  = FALSE;
167
168     t_forcetable        *pairsTable = nullptr; /* for 1-4 interactions, [pairs] and [pairs_nb] */
169
170     /* Free energy */
171     int      efep          = 0;
172     real     sc_alphavdw   = 0;
173     real     sc_alphacoul  = 0;
174     int      sc_power      = 0;
175     real     sc_r_power    = 0;
176     real     sc_sigma6_def = 0;
177     real     sc_sigma6_min = 0;
178
179     /* Information about atom properties for the molecule blocks in the system */
180     struct cginfo_mb_t *cginfo_mb                    = nullptr;
181     /* Information about atom properties for local and non-local atoms */
182     std::vector<int>    cginfo;
183
184     rvec               *shift_vec                    = nullptr;
185
186     int                 cutoff_scheme = 0;     /* group- or Verlet-style cutoff */
187     gmx_bool            bNonbonded    = FALSE; /* true if nonbonded calculations are *not* turned off */
188
189     /* The Nbnxm Verlet non-bonded machinery */
190     std::unique_ptr<nonbonded_verlet_t> nbv;
191
192     /* The wall tables (if used) */
193     int                    nwall    = 0;
194     t_forcetable        ***wall_tab = nullptr;
195
196     /* The number of atoms participating in do_force_lowlevel */
197     int natoms_force = 0;
198     /* The number of atoms participating in force and constraints */
199     int natoms_force_constr = 0;
200     /* The allocation size of vectors of size natoms_force */
201     int nalloc_force = 0;
202
203     /* Forces that should not enter into the coord x force virial summation:
204      * PPPM/PME/Ewald/posres/ForceProviders
205      */
206     /* True when we have contributions that are directly added to the virial */
207     bool                   haveDirectVirialContributions = false;
208     /* Force buffer for force computation with direct virial contributions */
209     std::vector<gmx::RVec> forceBufferForDirectVirialContributions;
210
211     /* Data for PPPM/PME/Ewald */
212     struct gmx_pme_t *pmedata                = nullptr;
213     int               ljpme_combination_rule = 0;
214
215     /* PME/Ewald stuff */
216     struct gmx_ewald_tab_t *ewald_table = nullptr;
217
218     /* Shift force array for computing the virial, size SHIFTS */
219     std::vector<gmx::RVec> shiftForces;
220
221     /* Non bonded Parameter lists */
222     int      ntype        = 0; /* Number of atom types */
223     gmx_bool bBHAM        = FALSE;
224     real    *nbfp         = nullptr;
225     real    *ljpme_c6grid = nullptr; /* C6-values used on grid in LJPME */
226
227     /* Energy group pair flags */
228     int *egp_flags = nullptr;
229
230     /* Shell molecular dynamics flexible constraints */
231     real fc_stepsize = 0;
232
233     /* If > 0 signals Test Particle Insertion,
234      * the value is the number of atoms of the molecule to insert
235      * Only the energy difference due to the addition of the last molecule
236      * should be calculated.
237      */
238     int n_tpi = 0;
239
240     /* QMMM stuff */
241     gmx_bool          bQMMM = FALSE;
242     struct t_QMMMrec *qr    = nullptr;
243
244     /* QM-MM neighborlists */
245     struct t_nblist        *QMMMlist = nullptr;
246
247     /* Limit for printing large forces, negative is don't print */
248     real print_force = 0;
249
250     /* User determined parameters, copied from the inputrec */
251     int  userint1  = 0;
252     int  userint2  = 0;
253     int  userint3  = 0;
254     int  userint4  = 0;
255     real userreal1 = 0;
256     real userreal2 = 0;
257     real userreal3 = 0;
258     real userreal4 = 0;
259
260     /* Pointer to struct for managing threading of bonded force calculation */
261     struct bonded_threading_t *bondedThreading = nullptr;
262
263     /* TODO: Replace the pointer by an object once we got rid of C */
264     gmx::GpuBonded *gpuBonded = nullptr;
265
266     /* Ewald correction thread local virial and energy data */
267     int                         nthread_ewc = 0;
268     struct ewald_corr_thread_t *ewc_t       = nullptr;
269
270     gmx::ForceProviders        *forceProviders = nullptr;
271 };
272
273 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
274  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
275  * in the code, but beware if you are using these macros externally.
276  */
277 #define C6(nbfp, ntp, ai, aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
278 #define C12(nbfp, ntp, ai, aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
279 #define BHAMC(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
280 #define BHAMA(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
281 #define BHAMB(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
282
283 #endif