Merge branch 'origin/release-2021' into merge-2021-into-master
[alexxy/gromacs.git] / src / gromacs / nbnxm / nbnxm.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 // FIXME: remove the "__" prefix in front of the group def when we move the
38 //        nonbonded code into separate dir.
39
40 /*! \libinternal \defgroup __module_nbnxm Short-range non-bonded interaction module
41  * \ingroup group_mdrun
42  *
43  * \brief Computes forces and energies for short-range pair-interactions
44  * based on the Verlet algorithm. The algorithm uses pair-lists generated
45  * at fixed intervals as well as various flavors of pair interaction kernels
46  * implemented for a wide range of CPU and GPU architectures.
47  *
48  * The module includes support for flavors of Coulomb and Lennard-Jones interaction
49  * treatment implemented for a large range of SIMD instruction sets for CPU
50  * architectures as well as in CUDA and OpenCL for GPU architectures.
51  * Additionally there is a reference CPU non-SIMD and a reference CPU
52  * for GPU pair-list setup interaction kernel.
53  *
54  * The implementation of the kernels is based on the cluster non-bonded algorithm
55  * which in the code is referred to as the NxM algorithms ("nbnxm_" prefix);
56  * for details of the algorithm see DOI:10.1016/j.cpc.2013.06.003.
57  *
58  * Algorithmically, the non-bonded computation has two different modes:
59  * A "classical" mode: generate a list every nstlist steps containing at least
60  * all atom pairs up to a distance of rlistOuter and compute pair interactions
61  * for all pairs that are within the interaction cut-off.
62  * A "dynamic pruning" mode: generate an "outer-list" up to cut-off rlistOuter
63  * every nstlist steps and prune the outer-list using a cut-off of rlistInner
64  * every nstlistPrune steps to obtain a, smaller, "inner-list". This
65  * results in fewer interaction computations and allows for a larger nstlist.
66  * On a GPU, this dynamic pruning is performed in a rolling fashion, pruning
67  * only a sub-part of the list each (second) step. This way it can often
68  * overlap with integration and constraints on the CPU.
69  * Currently a simple heuristic determines which mode will be used.
70  *
71  * TODO: add a summary list and brief descriptions of the different submodules:
72  * search, CPU kernels, GPU glue code + kernels.
73  *
74  * \author Berk Hess <hess@kth.se>
75  * \author Szilárd Páll <pall.szilard@gmail.com>
76  * \author Mark Abraham <mark.j.abraham@gmail.com>
77  * \author Anca Hamuraru <anca@streamcomputing.eu>
78  * \author Teemu Virolainen <teemu@streamcomputing.eu>
79  * \author Dimitrios Karkoulis <dimitris.karkoulis@gmail.com>
80  *
81  * TODO: add more authors!
82  */
83
84 /*! \libinternal
85  * \defgroup module_nbnxm Non-bonded pair interactions
86  * \ingroup group_mdrun
87  * \brief
88  * Implements non-bonded pair interaction functionality for NxM atom clusters.
89  *
90  * This module provides methods to, very efficiently, compute non-bonded
91  * pair interactions on CPUs as well as accelerators. It also provides
92  * a method to construct the NxM atom-cluster pair-list required for
93  * computing these non-bonded iteractions.
94  */
95
96 /*! \libinternal \file
97  *
98  * \brief This file contains the public interface of the nbnxm module
99  * that implements the NxM atom cluster non-bonded algorithm to efficiently
100  * compute pair forces.
101  *
102  *
103  * \author Berk Hess <hess@kth.se>
104  * \author Szilárd Páll <pall.szilard@gmail.com>
105  *
106  * \inlibraryapi
107  * \ingroup module_nbnxm
108  */
109
110
111 #ifndef GMX_NBNXM_NBNXM_H
112 #define GMX_NBNXM_NBNXM_H
113
114 #include <memory>
115
116 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
117 #include "gromacs/math/vectypes.h"
118 #include "gromacs/mdtypes/locality.h"
119 #include "gromacs/utility/arrayref.h"
120 #include "gromacs/utility/enumerationhelpers.h"
121 #include "gromacs/utility/real.h"
122
123 struct DeviceInformation;
124 struct gmx_domdec_zones_t;
125 struct gmx_enerdata_t;
126 struct gmx_hw_info_t;
127 struct gmx_mtop_t;
128 struct NbnxmGpu;
129 struct gmx_wallcycle;
130 struct interaction_const_t;
131 struct nbnxn_atomdata_t;
132 struct nonbonded_verlet_t;
133 class PairSearch;
134 class PairlistSets;
135 struct t_commrec;
136 struct t_lambda;
137 struct t_mdatoms;
138 struct t_nrnb;
139 struct t_forcerec;
140 struct t_inputrec;
141
142 class GpuEventSynchronizer;
143
144 namespace gmx
145 {
146 class DeviceStreamManager;
147 class ForceWithShiftForces;
148 class GpuBonded;
149 template<typename>
150 class ListOfLists;
151 class MDLogger;
152 template<typename>
153 class Range;
154 class StepWorkload;
155 class UpdateGroupsCog;
156 } // namespace gmx
157
158 //! Namespace for non-bonded kernels
159 namespace Nbnxm
160 {
161 enum class KernelType;
162
163 /*! \brief Nbnxm electrostatic GPU kernel flavors.
164  *
165  *  Types of electrostatics implementations available in the GPU non-bonded
166  *  force kernels. These represent both the electrostatics types implemented
167  *  by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
168  *  enums.h) as well as encode implementation details analytical/tabulated
169  *  and single or twin cut-off (for Ewald kernels).
170  *  Note that the cut-off and RF kernels have only analytical flavor and unlike
171  *  in the CPU kernels, the tabulated kernels are ATM Ewald-only.
172  *
173  *  The row-order of pointers to different electrostatic kernels defined in
174  *  nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
175  *  should match the order of enumerated types below.
176  */
177 enum class ElecType : int
178 {
179     Cut,          //!< Plain cut-off
180     RF,           //!< Reaction field
181     EwaldTab,     //!< Tabulated Ewald with single cut-off
182     EwaldTabTwin, //!< Tabulated Ewald with twin cut-off
183     EwaldAna,     //!< Analytical Ewald with single cut-off
184     EwaldAnaTwin, //!< Analytical Ewald with twin cut-off
185     Count         //!< Number of valid values
186 };
187
188 //! Number of possible \ref ElecType values.
189 constexpr int c_numElecTypes = static_cast<int>(ElecType::Count);
190
191 /*! \brief Nbnxm VdW GPU kernel flavors.
192  *
193  * The enumerates values correspond to the LJ implementations in the GPU non-bonded
194  * kernels.
195  *
196  * The column-order of pointers to different electrostatic kernels defined in
197  * nbnxn_cuda_ocl.cpp/.cu by the nb_*_kfunc_ptr function pointer table
198  * should match the order of enumerated types below.
199  */
200 enum class VdwType : int
201 {
202     Cut,         //!< Plain cut-off
203     CutCombGeom, //!< Cut-off with geometric combination rules
204     CutCombLB,   //!< Cut-off with Lorentz-Berthelot combination rules
205     FSwitch,     //!< Smooth force switch
206     PSwitch,     //!< Smooth potential switch
207     EwaldGeom,   //!< Ewald with geometric combination rules
208     EwaldLB,     //!< Ewald with Lorentz-Berthelot combination rules
209     Count        //!< Number of valid values
210 };
211
212 //! Number of possible \ref VdwType values.
213 constexpr int c_numVdwTypes = static_cast<int>(VdwType::Count);
214
215 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
216 enum class KernelType : int
217 {
218     NotSet = 0,
219     Cpu4x4_PlainC,
220     Cpu4xN_Simd_4xN,
221     Cpu4xN_Simd_2xNN,
222     Gpu8x8x8,
223     Cpu8x8x8_PlainC,
224     Count
225 };
226
227 /*! \brief Ewald exclusion types */
228 enum class EwaldExclusionType : int
229 {
230     NotSet = 0,
231     Table,
232     Analytical,
233     DecidedByGpuModule
234 };
235
236 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
237 struct KernelSetup
238 {
239     //! The non-bonded type, also affects the pairlist construction kernel
240     KernelType kernelType = KernelType::NotSet;
241     //! Ewald exclusion computation handling type, currently only used for CPU
242     EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
243 };
244
245 /*! \brief Return a string identifying the kernel type.
246  *
247  * \param [in] kernelType   nonbonded kernel type, takes values from the nbnxn_kernel_type enum
248  * \returns                 a string identifying the kernel corresponding to the type passed as argument
249  */
250 const char* lookup_kernel_name(Nbnxm::KernelType kernelType);
251
252 } // namespace Nbnxm
253
254 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
255 enum
256 {
257     enbvClearFNo,
258     enbvClearFYes
259 };
260
261 /*! \libinternal
262  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
263 struct nonbonded_verlet_t
264 {
265 public:
266     //! Constructs an object from its components
267     nonbonded_verlet_t(std::unique_ptr<PairlistSets>     pairlistSets,
268                        std::unique_ptr<PairSearch>       pairSearch,
269                        std::unique_ptr<nbnxn_atomdata_t> nbat,
270                        const Nbnxm::KernelSetup&         kernelSetup,
271                        NbnxmGpu*                         gpu_nbv,
272                        gmx_wallcycle*                    wcycle);
273
274     ~nonbonded_verlet_t();
275
276     //! Returns whether a GPU is use for the non-bonded calculations
277     bool useGpu() const { return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8; }
278
279     //! Returns whether a GPU is emulated for the non-bonded calculations
280     bool emulateGpu() const
281     {
282         return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
283     }
284
285     //! Return whether the pairlist is of simple, CPU type
286     bool pairlistIsSimple() const { return !useGpu() && !emulateGpu(); }
287
288
289     //! Returns the order of the local atoms on the grid
290     gmx::ArrayRef<const int> getLocalAtomOrder() const;
291
292     //! Sets the order of the local atoms to the order grid atom ordering
293     void setLocalAtomOrder();
294
295     //! Returns the index position of the atoms on the search grid
296     gmx::ArrayRef<const int> getGridIndices() const;
297
298     /*! \brief Constructs the pairlist for the given locality
299      *
300      * When there are no non-self exclusions, \p exclusions can be empty.
301      * Otherwise the number of lists in \p exclusions should match the number
302      * of atoms when not using DD, or the total number of atoms in the i-zones
303      * when using DD.
304      *
305      * \param[in] iLocality   The interaction locality: local or non-local
306      * \param[in] exclusions  Lists of exclusions for every atom.
307      * \param[in] step        Used to set the list creation step
308      * \param[in,out] nrnb    Flop accounting struct, can be nullptr
309      */
310     void constructPairlist(gmx::InteractionLocality     iLocality,
311                            const gmx::ListOfLists<int>& exclusions,
312                            int64_t                      step,
313                            t_nrnb*                      nrnb);
314
315     //! Updates all the atom properties in Nbnxm
316     void setAtomProperties(gmx::ArrayRef<const int>  atomTypes,
317                            gmx::ArrayRef<const real> atomCharges,
318                            gmx::ArrayRef<const int>  atomInfo);
319
320     /*!\brief Convert the coordinates to NBNXM format for the given locality.
321      *
322      * The API function for the transformation of the coordinates from one layout to another.
323      *
324      * \param[in] locality     Whether coordinates for local or non-local atoms should be
325      * transformed. \param[in] fillLocal    If the coordinates for filler particles should be
326      * zeroed. \param[in] coordinates  Coordinates in plain rvec format to be transformed.
327      */
328     void convertCoordinates(gmx::AtomLocality locality, bool fillLocal, gmx::ArrayRef<const gmx::RVec> coordinates);
329
330     /*!\brief Convert the coordinates to NBNXM format on the GPU for the given locality
331      *
332      * The API function for the transformation of the coordinates from one layout to another in the GPU memory.
333      *
334      * \param[in] locality        Whether coordinates for local or non-local atoms should be transformed.
335      * \param[in] fillLocal       If the coordinates for filler particles should be zeroed.
336      * \param[in] d_x             GPU coordinates buffer in plain rvec format to be transformed.
337      * \param[in] xReadyOnDevice  Event synchronizer indicating that the coordinates are ready in the device memory.
338      */
339     void convertCoordinatesGpu(gmx::AtomLocality       locality,
340                                bool                    fillLocal,
341                                DeviceBuffer<gmx::RVec> d_x,
342                                GpuEventSynchronizer*   xReadyOnDevice);
343
344     //! Init for GPU version of setup coordinates in Nbnxm
345     void atomdata_init_copy_x_to_nbat_x_gpu();
346
347     //! Sync the nonlocal GPU stream with dependent tasks in the local queue.
348     void insertNonlocalGpuDependency(gmx::InteractionLocality interactionLocality);
349
350     //! Returns a reference to the pairlist sets
351     const PairlistSets& pairlistSets() const { return *pairlistSets_; }
352
353     //! Returns whether step is a dynamic list pruning step, for CPU lists
354     bool isDynamicPruningStepCpu(int64_t step) const;
355
356     //! Returns whether step is a dynamic list pruning step, for GPU lists
357     bool isDynamicPruningStepGpu(int64_t step) const;
358
359     //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
360     void dispatchPruneKernelCpu(gmx::InteractionLocality iLocality, const rvec* shift_vec);
361
362     //! Dispatches the dynamic pruning kernel for GPU lists
363     void dispatchPruneKernelGpu(int64_t step);
364
365     //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
366     void dispatchNonbondedKernel(gmx::InteractionLocality   iLocality,
367                                  const interaction_const_t& ic,
368                                  const gmx::StepWorkload&   stepWork,
369                                  int                        clearF,
370                                  const t_forcerec&          fr,
371                                  gmx_enerdata_t*            enerd,
372                                  t_nrnb*                    nrnb);
373
374     //! Executes the non-bonded free-energy kernel, always runs on the CPU
375     void dispatchFreeEnergyKernel(gmx::InteractionLocality   iLocality,
376                                   const t_forcerec*          fr,
377                                   rvec                       x[],
378                                   gmx::ForceWithShiftForces* forceWithShiftForces,
379                                   const t_mdatoms&           mdatoms,
380                                   t_lambda*                  fepvals,
381                                   gmx::ArrayRef<const real>  lambda,
382                                   gmx_enerdata_t*            enerd,
383                                   const gmx::StepWorkload&   stepWork,
384                                   t_nrnb*                    nrnb);
385
386     /*! \brief Add the forces stored in nbat to f, zeros the forces in nbat
387      * \param [in] locality         Local or non-local
388      * \param [inout] force         Force to be added to
389      */
390     void atomdata_add_nbat_f_to_f(gmx::AtomLocality locality, gmx::ArrayRef<gmx::RVec> force);
391
392     /*! \brief Get the number of atoms for a given locality
393      *
394      * \param [in] locality   Local or non-local
395      * \returns               The number of atoms for given locality
396      */
397     int getNumAtoms(gmx::AtomLocality locality);
398
399     /*! \brief Get the pointer to the GPU nonbonded force buffer
400      *
401      * \returns A pointer to the force buffer in GPU memory
402      */
403     void* getGpuForces();
404
405     //! Return the kernel setup
406     const Nbnxm::KernelSetup& kernelSetup() const { return kernelSetup_; }
407
408     //! Returns the outer radius for the pair list
409     real pairlistInnerRadius() const;
410
411     //! Returns the outer radius for the pair list
412     real pairlistOuterRadius() const;
413
414     //! Changes the pair-list outer and inner radius
415     void changePairlistRadii(real rlistOuter, real rlistInner);
416
417     //! Set up internal flags that indicate what type of short-range work there is.
418     void setupGpuShortRangeWork(const gmx::GpuBonded* gpuBonded, gmx::InteractionLocality iLocality);
419
420     // TODO: Make all data members private
421 public:
422     //! All data related to the pair lists
423     std::unique_ptr<PairlistSets> pairlistSets_;
424     //! Working data for constructing the pairlists
425     std::unique_ptr<PairSearch> pairSearch_;
426     //! Atom data
427     std::unique_ptr<nbnxn_atomdata_t> nbat;
428
429 private:
430     //! The non-bonded setup, also affects the pairlist construction kernel
431     Nbnxm::KernelSetup kernelSetup_;
432     //! \brief Pointer to wallcycle structure.
433     gmx_wallcycle* wcycle_;
434
435 public:
436     //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
437     NbnxmGpu* gpu_nbv;
438 };
439
440 namespace Nbnxm
441 {
442
443 /*! \brief Creates an Nbnxm object */
444 std::unique_ptr<nonbonded_verlet_t> init_nb_verlet(const gmx::MDLogger& mdlog,
445                                                    const t_inputrec*    ir,
446                                                    const t_forcerec*    fr,
447                                                    const t_commrec*     cr,
448                                                    const gmx_hw_info_t& hardwareInfo,
449                                                    bool                 useGpuForNonbonded,
450                                                    const gmx::DeviceStreamManager* deviceStreamManager,
451                                                    const gmx_mtop_t*               mtop,
452                                                    matrix                          box,
453                                                    gmx_wallcycle*                  wcycle);
454
455 } // namespace Nbnxm
456
457 /*! \brief Put the atoms on the pair search grid.
458  *
459  * Only atoms with indices wihtin \p atomRange in x are put on the grid.
460  * When \p updateGroupsCog != nullptr, atoms are put on the grid
461  * based on the center of geometry of the group they belong to.
462  * Atoms or COGs of groups should be within the bounding box provided,
463  * this is checked in debug builds when not using update groups.
464  * The atom density is used to determine the grid size when \p gridIndex = 0.
465  * When \p atomDensity <= 0, the density is determined from atomEnd-atomStart
466  * and the bounding box corners.
467  * With domain decomposition, part of the atoms might have migrated,
468  * but have not been removed yet. This count is given by \p numAtomsMoved.
469  * When \p move[i] < 0 particle i has migrated and will not be put on the grid.
470  *
471  * \param[in,out] nb_verlet    The non-bonded object
472  * \param[in]     box          Box used for periodic distance calculations
473  * \param[in]     gridIndex    The index of the grid to spread to, always 0 except with test particle insertion
474  * \param[in]     lowerCorner  Atom groups to be gridded should have coordinates >= this corner
475  * \param[in]     upperCorner  Atom groups to be gridded should have coordinates <= this corner
476  * \param[in]     updateGroupsCog  Centers of geometry for update groups, pass nullptr when not using update groups
477  * \param[in]     atomRange    Range of atoms to grid
478  * \param[in]     atomDensity  An estimate of the atom density, used for peformance optimization and only with \p gridIndex = 0
479  * \param[in]     atomInfo     Atom information flags
480  * \param[in]     x            Coordinates for atoms to grid
481  * \param[in]     numAtomsMoved  The number of atoms that will move to another domain, pass 0 without DD
482  * \param[in]     move         Move flags for atoms, pass nullptr without DD
483  */
484 void nbnxn_put_on_grid(nonbonded_verlet_t*            nb_verlet,
485                        const matrix                   box,
486                        int                            gridIndex,
487                        const rvec                     lowerCorner,
488                        const rvec                     upperCorner,
489                        const gmx::UpdateGroupsCog*    updateGroupsCog,
490                        gmx::Range<int>                atomRange,
491                        real                           atomDensity,
492                        gmx::ArrayRef<const int>       atomInfo,
493                        gmx::ArrayRef<const gmx::RVec> x,
494                        int                            numAtomsMoved,
495                        const int*                     move);
496
497 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
498  *
499  * with domain decomposition. Should be called after calling
500  * nbnxn_search_put_on_grid for the local atoms / home zone.
501  */
502 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t*              nb_verlet,
503                                 const struct gmx_domdec_zones_t* zones,
504                                 gmx::ArrayRef<const int>         atomInfo,
505                                 gmx::ArrayRef<const gmx::RVec>   x);
506
507 #endif // GMX_NBNXN_NBNXM_H