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