f6404b48481645ea10aff18de59eb4cb27ec991b
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source 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 #ifndef GMX_MDTYPES_TYPES_FORCEREC_H
39 #define GMX_MDTYPES_TYPES_FORCEREC_H
40
41 #include <array>
42 #include <memory>
43 #include <vector>
44
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/real.h"
50
51 #include "locality.h"
52
53 /* Abstract type for PME that is defined only in the routine that use them. */
54 struct gmx_pme_t;
55 struct nonbonded_verlet_t;
56 struct bonded_threading_t;
57 class DeviceContext;
58 class DispersionCorrection;
59 class ListedForces;
60 struct t_fcdata;
61 struct t_forcetable;
62 struct interaction_const_t;
63
64 namespace gmx
65 {
66 class DeviceStreamManager;
67 class GpuBonded;
68 class GpuForceReduction;
69 class ForceProviders;
70 class StatePropagatorDataGpu;
71 class PmePpCommGpu;
72 class WholeMoleculeTransform;
73 } // namespace gmx
74
75 /* macros for the cginfo data in forcerec
76  *
77  * Since the tpx format support max 256 energy groups, we do the same here.
78  * Note that we thus have bits 8-14 still unused.
79  *
80  * The maximum cg size in cginfo is 63
81  * because we only have space for 6 bits in cginfo,
82  * this cg size entry is actually only read with domain decomposition.
83  */
84 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
85 #define GET_CGINFO_GID(cgi) ((cgi)&255)
86 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1 << 15))
87 #define GET_CGINFO_FEP(cgi) ((cgi) & (1 << 15))
88 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1 << 17))
89 #define GET_CGINFO_EXCL_INTER(cgi) ((cgi) & (1 << 17))
90 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1 << 20))
91 #define GET_CGINFO_CONSTR(cgi) ((cgi) & (1 << 20))
92 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1 << 21))
93 #define GET_CGINFO_SETTLE(cgi) ((cgi) & (1 << 21))
94 /* This bit is only used with bBondComm in the domain decomposition */
95 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1 << 22))
96 #define GET_CGINFO_BOND_INTER(cgi) ((cgi) & (1 << 22))
97 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1 << 23))
98 #define GET_CGINFO_HAS_VDW(cgi) ((cgi) & (1 << 23))
99 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1 << 24))
100 #define GET_CGINFO_HAS_Q(cgi) ((cgi) & (1 << 24))
101
102
103 /* Value to be used in mdrun for an infinite cut-off.
104  * Since we need to compare with the cut-off squared,
105  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
106  */
107 #define GMX_CUTOFF_INF 1E+18
108
109 /* enums for the neighborlist type */
110 enum
111 {
112     enbvdwNONE,
113     enbvdwLJ,
114     enbvdwBHAM,
115     enbvdwTAB,
116     enbvdwNR
117 };
118
119 struct cginfo_mb_t
120 {
121     int              cg_start = 0;
122     int              cg_end   = 0;
123     int              cg_mod   = 0;
124     std::vector<int> cginfo;
125 };
126
127
128 /* Forward declaration of type for managing Ewald tables */
129 struct gmx_ewald_tab_t;
130
131 struct ewald_corr_thread_t;
132
133 /*! \brief Helper force buffers for ForceOutputs
134  *
135  * This class stores intermediate force buffers that are used
136  * internally in the force calculation and which are reduced into
137  * the output force buffer passed to the force calculation.
138  */
139 class ForceHelperBuffers
140 {
141 public:
142     /*! \brief Constructs helper buffers
143      *
144      * When the forces that will be accumulated with help of these buffers
145      * have direct virial contributions, set the parameter to true, so
146      * an extra force buffer is available for these forces to enable
147      * correct virial computation.
148      */
149     ForceHelperBuffers(bool haveDirectVirialContributions);
150
151     //! Returns whether we have a direct virial contribution force buffer
152     bool haveDirectVirialContributions() const { return haveDirectVirialContributions_; }
153
154     //! Returns the buffer for direct virial contributions
155     gmx::ArrayRef<gmx::RVec> forceBufferForDirectVirialContributions()
156     {
157         GMX_ASSERT(haveDirectVirialContributions_, "Buffer can only be requested when present");
158         return forceBufferForDirectVirialContributions_;
159     }
160
161     //! Returns the buffer for shift forces, size SHIFTS
162     gmx::ArrayRef<gmx::RVec> shiftForces() { return shiftForces_; }
163
164     //! Resizes the direct virial contribution buffer, when present
165     void resize(int numAtoms);
166
167 private:
168     //! True when we have contributions that are directly added to the virial
169     bool haveDirectVirialContributions_ = false;
170     //! Force buffer for force computation with direct virial contributions
171     std::vector<gmx::RVec> forceBufferForDirectVirialContributions_;
172     //! Shift force array for computing the virial, size SHIFTS
173     std::vector<gmx::RVec> shiftForces_;
174 };
175 // NOLINTNEXTLINE (clang-analyzer-optin.performance.Padding)
176 struct t_forcerec
177 {
178     // Declare an explicit constructor and destructor, so they can be
179     // implemented in a single source file, so that not every source
180     // file that includes this one needs to understand how to find the
181     // destructors of the objects pointed to by unique_ptr members.
182     t_forcerec();
183     ~t_forcerec();
184
185     std::unique_ptr<interaction_const_t> ic;
186
187     /* PBC stuff */
188     PbcType pbcType = PbcType::Xyz;
189     //! Tells whether atoms inside a molecule can be in different periodic images,
190     //  i.e. whether we need to take into account PBC when computing distances inside molecules.
191     //  This determines whether PBC must be considered for e.g. bonded interactions.
192     gmx_bool        bMolPBC     = FALSE;
193     RefCoordScaling rc_scaling  = RefCoordScaling::No;
194     rvec            posres_com  = { 0 };
195     rvec            posres_comB = { 0 };
196
197     gmx_bool use_simd_kernels = FALSE;
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     NbkernelElecType     nbkernel_elec_interaction = NbkernelElecType::None;
209     NbkernelVdwType      nbkernel_vdw_interaction  = NbkernelVdwType::None;
210     InteractionModifiers nbkernel_elec_modifier    = InteractionModifiers::None;
211     InteractionModifiers nbkernel_vdw_modifier     = InteractionModifiers::None;
212
213     /* Cut-Off stuff.
214      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
215      */
216     real rlist = 0;
217
218     /* Charge sum for topology A/B ([0]/[1]) for Ewald corrections */
219     double qsum[2]  = { 0 };
220     double q2sum[2] = { 0 };
221     double c6sum[2] = { 0 };
222
223     /* Dispersion correction stuff */
224     std::unique_ptr<DispersionCorrection> dispersionCorrection;
225
226     /* Fudge factors */
227     real fudgeQQ = 0;
228
229     /* Table stuff */
230     gmx_bool bcoultab = FALSE;
231     gmx_bool bvdwtab  = FALSE;
232
233     std::unique_ptr<t_forcetable> pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
234
235     /* Free energy */
236     FreeEnergyPerturbationType efep = FreeEnergyPerturbationType::No;
237
238     /* Information about atom properties for the molecule blocks in the system */
239     std::vector<cginfo_mb_t> cginfo_mb;
240     /* Information about atom properties for local and non-local atoms */
241     std::vector<int> cginfo;
242
243     rvec* shift_vec = nullptr;
244
245     std::unique_ptr<gmx::WholeMoleculeTransform> wholeMoleculeTransform;
246
247     /* The Nbnxm Verlet non-bonded machinery */
248     std::unique_ptr<nonbonded_verlet_t> nbv;
249
250     /* The wall tables (if used) */
251     int                                                     nwall = 0;
252     std::vector<std::vector<std::unique_ptr<t_forcetable>>> wall_tab;
253
254     /* The number of atoms participating in do_force_lowlevel */
255     int natoms_force = 0;
256     /* The number of atoms participating in force calculation and constraints */
257     int natoms_force_constr = 0;
258
259     /* List of helper buffers for ForceOutputs, one for each time step with MTS */
260     std::vector<ForceHelperBuffers> forceHelperBuffers;
261
262     /* Data for PPPM/PME/Ewald */
263     struct gmx_pme_t* pmedata                = nullptr;
264     LongRangeVdW      ljpme_combination_rule = LongRangeVdW::Geom;
265
266     /* PME/Ewald stuff */
267     std::unique_ptr<gmx_ewald_tab_t> ewald_table;
268
269     /* Non bonded Parameter lists */
270     int               ntype = 0; /* Number of atom types */
271     gmx_bool          bBHAM = FALSE;
272     std::vector<real> nbfp;
273     std::vector<real> ljpme_c6grid; /* C6-values used on grid in LJPME */
274
275     /* Energy group pair flags */
276     int* egp_flags = nullptr;
277
278     /* Shell molecular dynamics flexible constraints */
279     real fc_stepsize = 0;
280
281     /* If > 0 signals Test Particle Insertion,
282      * the value is the number of atoms of the molecule to insert
283      * Only the energy difference due to the addition of the last molecule
284      * should be calculated.
285      */
286     int n_tpi = 0;
287
288     /* Limit for printing large forces, negative is don't print */
289     real print_force = 0;
290
291     /* User determined parameters, copied from the inputrec */
292     int  userint1  = 0;
293     int  userint2  = 0;
294     int  userint3  = 0;
295     int  userint4  = 0;
296     real userreal1 = 0;
297     real userreal2 = 0;
298     real userreal3 = 0;
299     real userreal4 = 0;
300
301     /* Tells whether we use multiple time stepping, computing some forces less frequently */
302     bool useMts = false;
303
304     /* Data for special listed force calculations */
305     std::unique_ptr<t_fcdata> fcdata;
306
307     // The listed forces calculation data, 1 entry or multiple entries with multiple time stepping
308     std::vector<ListedForces> listedForces;
309
310     /* TODO: Replace the pointer by an object once we got rid of C */
311     gmx::GpuBonded* gpuBonded = nullptr;
312
313     /* Ewald correction thread local virial and energy data */
314     int                         nthread_ewc = 0;
315     struct ewald_corr_thread_t* ewc_t       = nullptr;
316
317     gmx::ForceProviders* forceProviders = nullptr;
318
319     // The stateGpu object is created in runner, forcerec just keeps the copy of the pointer.
320     // TODO: This is not supposed to be here. StatePropagatorDataGpu should be a part of
321     //       general StatePropagatorData object that is passed around
322     gmx::StatePropagatorDataGpu* stateGpu = nullptr;
323     // TODO: Should not be here. This is here only to pass the pointer around.
324     gmx::DeviceStreamManager* deviceStreamManager = nullptr;
325
326     //! GPU device context
327     DeviceContext* deviceContext = nullptr;
328
329     /* For PME-PP GPU communication */
330     std::unique_ptr<gmx::PmePpCommGpu> pmePpCommGpu;
331
332     /* For GPU force reduction (on both local and non-local atoms) */
333     gmx::EnumerationArray<gmx::AtomLocality, std::unique_ptr<gmx::GpuForceReduction>> gpuForceReduction;
334 };
335
336 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
337  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
338  * in the code, but beware if you are using these macros externally.
339  */
340 #define C6(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj))]
341 #define C12(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj)) + 1]
342 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj))]
343 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 1]
344 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 2]
345
346 #endif