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