Partial renaming of cginfo to atom info
[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/ishift.h"
48 #include "gromacs/pbcutil/pbc.h"
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/real.h"
51
52 #include "locality.h"
53
54 /* Abstract type for PME that is defined only in the routine that use them. */
55 struct gmx_pme_t;
56 struct nonbonded_verlet_t;
57 struct bonded_threading_t;
58 class DeviceContext;
59 class DispersionCorrection;
60 class ListedForces;
61 struct t_fcdata;
62 struct t_forcetable;
63 struct interaction_const_t;
64
65 namespace gmx
66 {
67 class DeviceStreamManager;
68 class ListedForcesGpu;
69 class GpuForceReduction;
70 class ForceProviders;
71 class StatePropagatorDataGpu;
72 class PmePpCommGpu;
73 class WholeMoleculeTransform;
74 } // namespace gmx
75
76 /* macros for the atom info data in forcerec
77  *
78  * Since the tpx format support max 256 energy groups, we do the same here.
79  * Note that we thus have bits 8-14 still unused.
80  */
81 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
82 #define GET_CGINFO_GID(cgi) ((cgi)&255)
83 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1 << 15))
84 #define GET_CGINFO_FEP(cgi) ((cgi) & (1 << 15))
85 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1 << 17))
86 #define GET_CGINFO_EXCL_INTER(cgi) ((cgi) & (1 << 17))
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
99
100 /* Value to be used in mdrun for an infinite cut-off.
101  * Since we need to compare with the cut-off squared,
102  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
103  */
104 #define GMX_CUTOFF_INF 1E+18
105 //! Check the cuttoff
106 real cutoff_inf(real cutoff);
107
108 /*! \brief Contains information about each atom in a molecule block of
109  * the global topology. */
110 struct AtomInfoWithinMoleculeBlock
111 {
112     //! Index within the system of the first atom in the molecule block
113     int indexOfFirstAtomInMoleculeBlock = 0;
114     //! Index within the system of the last atom in the molecule block
115     int indexOfLastAtomInMoleculeBlock = 0;
116     /*! \brief Atom info for each atom in the block.
117      *
118      * The typical case is that all atoms are identical for each
119      * molecule of the block, and if so this vector has size equal to
120      * the number of atoms in the molecule.
121      *
122      * An example of an atypical case is QM/MM, where multiple
123      * molecules might be present and different molecules have
124      * different atoms within any one QM group or region. Now there are
125      * multiple kinds of molecules with the same connectivity, so we simply
126      * write out the atom info for the entire molecule block. Then the
127      * size equals the product of the number of atoms in the
128      * molecule and the number of molecules.
129      *
130      * The vector needs to be indexed accordingly. */
131     std::vector<int> atomInfo;
132 };
133
134
135 /* Forward declaration of type for managing Ewald tables */
136 struct gmx_ewald_tab_t;
137
138 struct ewald_corr_thread_t;
139
140 /*! \brief Helper force buffers for ForceOutputs
141  *
142  * This class stores intermediate force buffers that are used
143  * internally in the force calculation and which are reduced into
144  * the output force buffer passed to the force calculation.
145  */
146 class ForceHelperBuffers
147 {
148 public:
149     /*! \brief Constructs helper buffers
150      *
151      * When the forces that will be accumulated with help of these buffers
152      * have direct virial contributions, set the parameter to true, so
153      * an extra force buffer is available for these forces to enable
154      * correct virial computation.
155      */
156     ForceHelperBuffers(bool haveDirectVirialContributions);
157
158     //! Returns whether we have a direct virial contribution force buffer
159     bool haveDirectVirialContributions() const { return haveDirectVirialContributions_; }
160
161     //! Returns the buffer for direct virial contributions
162     gmx::ArrayRef<gmx::RVec> forceBufferForDirectVirialContributions()
163     {
164         GMX_ASSERT(haveDirectVirialContributions_, "Buffer can only be requested when present");
165         return forceBufferForDirectVirialContributions_;
166     }
167
168     //! Returns the buffer for shift forces, size c_numShiftVectors
169     gmx::ArrayRef<gmx::RVec> shiftForces() { return shiftForces_; }
170
171     //! Resizes the direct virial contribution buffer, when present
172     void resize(int numAtoms);
173
174 private:
175     //! True when we have contributions that are directly added to the virial
176     bool haveDirectVirialContributions_ = false;
177     //! Force buffer for force computation with direct virial contributions
178     std::vector<gmx::RVec> forceBufferForDirectVirialContributions_;
179     //! Shift force array for computing the virial, size c_numShiftVectors
180     std::vector<gmx::RVec> shiftForces_;
181 };
182 // NOLINTNEXTLINE (clang-analyzer-optin.performance.Padding)
183 struct t_forcerec
184 {
185     // Declare an explicit constructor and destructor, so they can be
186     // implemented in a single source file, so that not every source
187     // file that includes this one needs to understand how to find the
188     // destructors of the objects pointed to by unique_ptr members.
189     t_forcerec();
190     ~t_forcerec();
191
192     std::unique_ptr<interaction_const_t> ic;
193
194     /* PBC stuff */
195     PbcType pbcType = PbcType::Xyz;
196     //! Tells whether atoms inside a molecule can be in different periodic images,
197     //  i.e. whether we need to take into account PBC when computing distances inside molecules.
198     //  This determines whether PBC must be considered for e.g. bonded interactions.
199     bool            bMolPBC     = false;
200     RefCoordScaling rc_scaling  = RefCoordScaling::No;
201     gmx::RVec       posres_com  = { 0, 0, 0 };
202     gmx::RVec       posres_comB = { 0, 0, 0 };
203
204     bool use_simd_kernels = false;
205
206     /* Interaction for calculated in kernels. In many cases this is similar to
207      * the electrostatics settings in the inputrecord, but the difference is that
208      * these variables always specify the actual interaction in the kernel - if
209      * we are tabulating reaction-field the inputrec will say reaction-field, but
210      * the kernel interaction will say cubic-spline-table. To be safe we also
211      * have a kernel-specific setting for the modifiers - if the interaction is
212      * tabulated we already included the inputrec modification there, so the kernel
213      * modification setting will say 'none' in that case.
214      */
215     NbkernelElecType     nbkernel_elec_interaction = NbkernelElecType::None;
216     NbkernelVdwType      nbkernel_vdw_interaction  = NbkernelVdwType::None;
217     InteractionModifiers nbkernel_elec_modifier    = InteractionModifiers::None;
218     InteractionModifiers nbkernel_vdw_modifier     = InteractionModifiers::None;
219
220     /* Cut-Off stuff.
221      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
222      */
223     real rlist = 0;
224
225     /* Charge sum for topology A/B ([0]/[1]) for Ewald corrections */
226     std::array<double, 2> qsum  = { 0 };
227     std::array<double, 2> q2sum = { 0 };
228     std::array<double, 2> c6sum = { 0 };
229
230     /* Dispersion correction stuff */
231     std::unique_ptr<DispersionCorrection> dispersionCorrection;
232
233     /* Fudge factors */
234     real fudgeQQ = 0;
235
236     std::unique_ptr<t_forcetable> pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
237
238     /* Free energy */
239     FreeEnergyPerturbationType efep = FreeEnergyPerturbationType::No;
240
241     /* Information about atom properties for the molecule blocks in the global topology */
242     std::vector<AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock;
243     /* Information about atom properties for local and non-local atoms */
244     std::vector<int> atomInfo;
245
246     std::vector<gmx::RVec> shift_vec;
247
248     std::unique_ptr<gmx::WholeMoleculeTransform> wholeMoleculeTransform;
249
250     /* The Nbnxm Verlet non-bonded machinery */
251     std::unique_ptr<nonbonded_verlet_t> nbv;
252
253     /* The wall tables (if used) */
254     int                                                     nwall = 0;
255     std::vector<std::vector<std::unique_ptr<t_forcetable>>> wall_tab;
256
257     /* The number of atoms participating in do_force_lowlevel */
258     int natoms_force = 0;
259     /* The number of atoms participating in force calculation and constraints */
260     int natoms_force_constr = 0;
261
262     /* List of helper buffers for ForceOutputs, one for each time step with MTS */
263     std::vector<ForceHelperBuffers> forceHelperBuffers;
264
265     /* Data for PPPM/PME/Ewald */
266     gmx_pme_t*   pmedata                = nullptr;
267     LongRangeVdW ljpme_combination_rule = LongRangeVdW::Geom;
268
269     /* PME/Ewald stuff */
270     std::unique_ptr<gmx_ewald_tab_t> ewald_table;
271
272     /* Non bonded Parameter lists */
273     int               ntype          = 0; /* Number of atom types */
274     bool              haveBuckingham = false;
275     std::vector<real> nbfp;
276     std::vector<real> ljpme_c6grid; /* C6-values used on grid in LJPME */
277
278     /* Energy group pair flags */
279     int* egp_flags = nullptr;
280
281     /* Shell molecular dynamics flexible constraints */
282     real fc_stepsize = 0;
283
284     /* If > 0 signals Test Particle Insertion,
285      * the value is the number of atoms of the molecule to insert
286      * Only the energy difference due to the addition of the last molecule
287      * should be calculated.
288      */
289     int n_tpi = 0;
290
291     /* Limit for printing large forces, negative is don't print */
292     real print_force = 0;
293
294     /* User determined parameters, copied from the inputrec */
295     int  userint1  = 0;
296     int  userint2  = 0;
297     int  userint3  = 0;
298     int  userint4  = 0;
299     real userreal1 = 0;
300     real userreal2 = 0;
301     real userreal3 = 0;
302     real userreal4 = 0;
303
304     /* Tells whether we use multiple time stepping, computing some forces less frequently */
305     bool useMts = false;
306
307     /* Data for special listed force calculations */
308     std::unique_ptr<t_fcdata> fcdata;
309
310     // The listed forces calculation data, 1 entry or multiple entries with multiple time stepping
311     std::vector<ListedForces> listedForces;
312
313     /* TODO: Replace the pointer by an object once we got rid of C */
314     gmx::ListedForcesGpu* listedForcesGpu = nullptr;
315
316     /* Ewald correction thread local virial and energy data */
317     int                              nthread_ewc = 0;
318     std::vector<ewald_corr_thread_t> ewc_t;
319
320     gmx::ForceProviders* forceProviders = nullptr;
321
322     // The stateGpu object is created in runner, forcerec just keeps the copy of the pointer.
323     // TODO: This is not supposed to be here. StatePropagatorDataGpu should be a part of
324     //       general StatePropagatorData object that is passed around
325     gmx::StatePropagatorDataGpu* stateGpu = nullptr;
326     // TODO: Should not be here. This is here only to pass the pointer around.
327     gmx::DeviceStreamManager* deviceStreamManager = nullptr;
328
329     //! GPU device context
330     DeviceContext* deviceContext = nullptr;
331
332     /* For PME-PP GPU communication */
333     std::unique_ptr<gmx::PmePpCommGpu> pmePpCommGpu;
334
335     /* For GPU force reduction (on both local and non-local atoms) */
336     gmx::EnumerationArray<gmx::AtomLocality, std::unique_ptr<gmx::GpuForceReduction>> gpuForceReduction;
337 };
338
339 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
340  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
341  * in the code, but beware if you are using these macros externally.
342  */
343 #define C6(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj))]
344 #define C12(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj)) + 1]
345 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj))]
346 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 1]
347 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 2]
348
349 #endif