Use gmx::Range in Nbnxm gridding functions
[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
84  * \defgroup module_nbnxm Non-bonded pair interactions
85  * \ingroup group_mdrun
86  * \brief
87  * Implements non-bonded pair interaction functionality for NxM atom clusters.
88  *
89  * This module provides methods to, very efficiently, compute non-bonded
90  * pair interactions on CPUs as well as accelerators. It also provides
91  * a method to construct the NxM atom-cluster pair-list required for
92  * computing these non-bonded iteractions.
93  */
94
95 /*! \libinternal \file
96  *
97  * \brief This file contains the public interface of the nbnxm module
98  * that implements the NxM atom cluster non-bonded algorithm to efficiently
99  * compute pair forces.
100  *
101  *
102  * \author Berk Hess <hess@kth.se>
103  * \author Szilárd Páll <pall.szilard@gmail.com>
104  *
105  * \inlibraryapi
106  * \ingroup module_nbnxm
107  */
108
109
110 #ifndef GMX_NBNXM_NBNXM_H
111 #define GMX_NBNXM_NBNXM_H
112
113 #include <memory>
114
115 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
116 #include "gromacs/math/vectypes.h"
117 #include "gromacs/utility/arrayref.h"
118 #include "gromacs/utility/enumerationhelpers.h"
119 #include "gromacs/utility/range.h"
120 #include "gromacs/utility/real.h"
121
122 #include "locality.h"
123
124 // TODO: Remove this include and the two nbnxm includes above
125 #include "nbnxm_gpu.h"
126
127 struct gmx_device_info_t;
128 struct gmx_domdec_zones_t;
129 struct gmx_enerdata_t;
130 struct gmx_hw_info_t;
131 struct gmx_mtop_t;
132 struct gmx_wallcycle;
133 struct interaction_const_t;
134 struct nonbonded_verlet_t;
135 class PairSearch;
136 class PairlistSets;
137 struct t_blocka;
138 struct t_commrec;
139 struct t_lambda;
140 struct t_mdatoms;
141 struct t_nrnb;
142 struct t_forcerec;
143 struct t_inputrec;
144
145 /*! \brief Switch for whether to use GPU for buffer ops*/
146 enum class BufferOpsUseGpu
147 {
148     True,
149     False
150 };
151
152 class GpuEventSynchronizer;
153
154 namespace gmx
155 {
156 class ForceWithShiftForces;
157 class MDLogger;
158 class UpdateGroupsCog;
159 }
160
161 namespace Nbnxm
162 {
163 enum class KernelType;
164 }
165
166 namespace Nbnxm
167 {
168
169 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
170 enum class KernelType : int
171 {
172     NotSet = 0,
173     Cpu4x4_PlainC,
174     Cpu4xN_Simd_4xN,
175     Cpu4xN_Simd_2xNN,
176     Gpu8x8x8,
177     Cpu8x8x8_PlainC,
178     Count
179 };
180
181 /*! \brief Ewald exclusion types */
182 enum class EwaldExclusionType : int
183 {
184     NotSet = 0,
185     Table,
186     Analytical,
187     DecidedByGpuModule
188 };
189
190 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
191 struct KernelSetup
192 {
193     //! The non-bonded type, also affects the pairlist construction kernel
194     KernelType         kernelType = KernelType::NotSet;
195     //! Ewald exclusion computation handling type, currently only used for CPU
196     EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
197 };
198
199 /*! \brief Return a string identifying the kernel type.
200  *
201  * \param [in] kernelType   nonbonded kernel type, takes values from the nbnxn_kernel_type enum
202  * \returns                 a string identifying the kernel corresponding to the type passed as argument
203  */
204 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
205
206 } // namespace Nbnxm
207
208 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
209 enum {
210     enbvClearFNo, enbvClearFYes
211 };
212
213 /*! \libinternal
214  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
215 struct nonbonded_verlet_t
216 {
217     public:
218         //! Constructs an object from its components
219         nonbonded_verlet_t(std::unique_ptr<PairlistSets>      pairlistSets,
220                            std::unique_ptr<PairSearch>        pairSearch,
221                            std::unique_ptr<nbnxn_atomdata_t>  nbat,
222                            const Nbnxm::KernelSetup          &kernelSetup,
223                            gmx_nbnxn_gpu_t                   *gpu_nbv,
224                            gmx_wallcycle                     *wcycle);
225
226         ~nonbonded_verlet_t();
227
228         //! Returns whether a GPU is use for the non-bonded calculations
229         bool useGpu() const
230         {
231             return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
232         }
233
234         //! Returns whether a GPU is emulated for the non-bonded calculations
235         bool emulateGpu() const
236         {
237             return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
238         }
239
240         //! Return whether the pairlist is of simple, CPU type
241         bool pairlistIsSimple() const
242         {
243             return !useGpu() && !emulateGpu();
244         }
245
246         //! Initialize the pair list sets, TODO this should be private
247         void initPairlistSets(bool haveMultipleDomains);
248
249         //! Returns the order of the local atoms on the grid
250         gmx::ArrayRef<const int> getLocalAtomOrder() const;
251
252         //! Sets the order of the local atoms to the order grid atom ordering
253         void setLocalAtomOrder();
254
255         //! Returns the index position of the atoms on the search grid
256         gmx::ArrayRef<const int> getGridIndices() const;
257
258         //! Constructs the pairlist for the given locality
259         void constructPairlist(Nbnxm::InteractionLocality  iLocality,
260                                const t_blocka             *excl,
261                                int64_t                     step,
262                                t_nrnb                     *nrnb);
263
264         //! Updates all the atom properties in Nbnxm
265         void setAtomProperties(const t_mdatoms          &mdatoms,
266                                gmx::ArrayRef<const int>  atomInfo);
267
268         /*!\brief Convert the coordinates to NBNXM format for the given locality.
269          *
270          * The API function for the transformation of the coordinates from one layout to another.
271          *
272          * \param[in] locality     Whether coordinates for local or non-local atoms should be transformed.
273          * \param[in] fillLocal    If the coordinates for filler particles should be zeroed.
274          * \param[in] coordinates  Coordinates in plain rvec format to be transformed.
275          */
276         void convertCoordinates(Nbnxm::AtomLocality             locality,
277                                 bool                            fillLocal,
278                                 gmx::ArrayRef<const gmx::RVec>  coordinates);
279
280         /*!\brief Convert the coordinates to NBNXM format on the GPU for the given locality
281          *
282          * The API function for the transformation of the coordinates from one layout to another in the GPU memory.
283          *
284          * \param[in] locality        Whether coordinates for local or non-local atoms should be transformed.
285          * \param[in] fillLocal       If the coordinates for filler particles should be zeroed.
286          * \param[in] d_x             GPU coordinates buffer in plain rvec format to be transformed.
287          * \param[in] xReadyOnDevice  Event synchronizer indicating that the coordinates are ready in the device memory.
288          */
289         void convertCoordinatesGpu(Nbnxm::AtomLocality              locality,
290                                    bool                             fillLocal,
291                                    DeviceBuffer<float>              d_x,
292                                    GpuEventSynchronizer            *xReadyOnDevice);
293
294         //! Init for GPU version of setup coordinates in Nbnxm
295         void atomdata_init_copy_x_to_nbat_x_gpu();
296
297         //! Sync the nonlocal GPU stream with dependent tasks in the local queue.
298         void insertNonlocalGpuDependency(Nbnxm::InteractionLocality interactionLocality);
299
300         //! Returns a reference to the pairlist sets
301         const PairlistSets &pairlistSets() const
302         {
303             return *pairlistSets_;
304         }
305
306         //! Returns whether step is a dynamic list pruning step, for CPU lists
307         bool isDynamicPruningStepCpu(int64_t step) const;
308
309         //! Returns whether step is a dynamic list pruning step, for GPU lists
310         bool isDynamicPruningStepGpu(int64_t step) const;
311
312         //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
313         void dispatchPruneKernelCpu(Nbnxm::InteractionLocality  iLocality,
314                                     const rvec                 *shift_vec);
315
316         //! Dispatches the dynamic pruning kernel for GPU lists
317         void dispatchPruneKernelGpu(int64_t step);
318
319         //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
320         void dispatchNonbondedKernel(Nbnxm::InteractionLocality  iLocality,
321                                      const interaction_const_t  &ic,
322                                      const gmx::StepWorkload    &stepWork,
323                                      int                         clearF,
324                                      const t_forcerec           &fr,
325                                      gmx_enerdata_t             *enerd,
326                                      t_nrnb                     *nrnb);
327
328         //! Executes the non-bonded free-energy kernel, always runs on the CPU
329         void dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocality,
330                                       const t_forcerec           *fr,
331                                       rvec                        x[],
332                                       gmx::ForceWithShiftForces  *forceWithShiftForces,
333                                       const t_mdatoms            &mdatoms,
334                                       t_lambda                   *fepvals,
335                                       real                       *lambda,
336                                       gmx_enerdata_t             *enerd,
337                                       const gmx::StepWorkload    &stepWork,
338                                       t_nrnb                     *nrnb);
339
340         /*! \brief Add the forces stored in nbat to f, zeros the forces in nbat
341          * \param [in] locality         Local or non-local
342          * \param [inout] force         Force to be added to
343          */
344         void atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality                 locality,
345                                       gmx::ArrayRef<gmx::RVec>            force);
346
347         /*! \brief Add the forces stored in nbat to total force using GPU buffer opse
348          *
349          * \param [in]     locality             Local or non-local
350          * \param [in,out] totalForcesDevice    Force to be added to
351          * \param [in]     forcesPmeDevice      Device buffer with PME forces
352          * \param[in]      dependencyList       List of synchronizers that represent the dependencies the reduction task needs to sync on.
353          * \param [in]     useGpuFPmeReduction  Whether PME forces should be added
354          * \param [in]     accumulateForce      If the total force buffer already contains data
355          */
356         void atomdata_add_nbat_f_to_f_gpu(Nbnxm::AtomLocality                         locality,
357                                           DeviceBuffer<float>                         totalForcesDevice,
358                                           void                                       *forcesPmeDevice,
359                                           gmx::ArrayRef<GpuEventSynchronizer* const>  dependencyList,
360                                           bool                                        useGpuFPmeReduction,
361                                           bool                                        accumulateForce);
362
363         /*! \brief Outer body of function to perform initialization for F buffer operations on GPU.
364          *
365          * \param localReductionDone     Pointer to an event synchronizer that marks the completion of the local f buffer ops kernel.
366          */
367         void atomdata_init_add_nbat_f_to_f_gpu(GpuEventSynchronizer* localReductionDone);
368
369         /*! \brief return pointer to GPU event recorded when coordinates have been copied to device */
370         void* get_x_on_device_event();
371
372         /*! \brief Wait for non-local copy of coordinate buffer from device to host */
373         void wait_nonlocal_x_copy_D2H_done();
374
375         /*! \brief return GPU pointer to f in rvec format */
376         void* get_gpu_frvec();
377
378         /*! \brief Ensure local stream waits for non-local stream */
379         void stream_local_wait_for_nonlocal();
380
381         //! Return the kernel setup
382         const Nbnxm::KernelSetup &kernelSetup() const
383         {
384             return kernelSetup_;
385         }
386
387         //! Returns the outer radius for the pair list
388         real pairlistInnerRadius() const;
389
390         //! Returns the outer radius for the pair list
391         real pairlistOuterRadius() const;
392
393         //! Changes the pair-list outer and inner radius
394         void changePairlistRadii(real rlistOuter,
395                                  real rlistInner);
396
397         //! Set up internal flags that indicate what type of short-range work there is.
398         void setupGpuShortRangeWork(const gmx::GpuBonded             *gpuBonded,
399                                     const Nbnxm::InteractionLocality  iLocality)
400         {
401             if (useGpu() && !emulateGpu())
402             {
403                 Nbnxm::setupGpuShortRangeWork(gpu_nbv, gpuBonded, iLocality);
404             }
405         }
406
407         //! Returns true if there is GPU short-range work for the given atom locality.
408         bool haveGpuShortRangeWork(const Nbnxm::AtomLocality aLocality)
409         {
410             return ((useGpu() && !emulateGpu()) &&
411                     Nbnxm::haveGpuShortRangeWork(gpu_nbv, aLocality));
412         }
413
414         // TODO: Make all data members private
415     public:
416         //! All data related to the pair lists
417         std::unique_ptr<PairlistSets>     pairlistSets_;
418         //! Working data for constructing the pairlists
419         std::unique_ptr<PairSearch>       pairSearch_;
420         //! Atom data
421         std::unique_ptr<nbnxn_atomdata_t> nbat;
422     private:
423         //! The non-bonded setup, also affects the pairlist construction kernel
424         Nbnxm::KernelSetup                kernelSetup_;
425         //! \brief Pointer to wallcycle structure.
426         gmx_wallcycle                    *wcycle_;
427     public:
428         //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
429         gmx_nbnxn_gpu_t                  *gpu_nbv;
430 };
431
432 namespace Nbnxm
433 {
434
435 /*! \brief Creates an Nbnxm object */
436 std::unique_ptr<nonbonded_verlet_t>
437 init_nb_verlet(const gmx::MDLogger     &mdlog,
438                gmx_bool                 bFEP_NonBonded,
439                const t_inputrec        *ir,
440                const t_forcerec        *fr,
441                const t_commrec         *cr,
442                const gmx_hw_info_t     &hardwareInfo,
443                const gmx_device_info_t *deviceInfo,
444                const gmx_mtop_t        *mtop,
445                matrix                   box,
446                gmx_wallcycle           *wcycle);
447
448 } // namespace Nbnxm
449
450 /*! \brief Put the atoms on the pair search grid.
451  *
452  * Only atoms with indices wihtin \p atomRange in x are put on the grid.
453  * When \p updateGroupsCog != nullptr, atoms are put on the grid
454  * based on the center of geometry of the group they belong to.
455  * Atoms or COGs of groups should be within the bounding box provided,
456  * this is checked in debug builds when not using update groups.
457  * The atom density is used to determine the grid size when \p gridIndex = 0.
458  * When \p atomDensity <= 0, the density is determined from atomEnd-atomStart
459  * and the bounding box corners.
460  * With domain decomposition, part of the atoms might have migrated,
461  * but have not been removed yet. This count is given by \p numAtomsMoved.
462  * When \p move[i] < 0 particle i has migrated and will not be put on the grid.
463  *
464  * \param[in,out] nb_verlet    The non-bonded object
465  * \param[in]     box          Box used for periodic distance calculations
466  * \param[in]     gridIndex    The index of the grid to spread to, always 0 except with test particle insertion
467  * \param[in]     lowerCorner  Atom groups to be gridded should have coordinates >= this corner
468  * \param[in]     upperCorner  Atom groups to be gridded should have coordinates <= this corner
469  * \param[in]     updateGroupsCog  Centers of geometry for update groups, pass nullptr when not using update groups
470  * \param[in]     atomRange    Range of atoms to grid
471  * \param[in]     atomDensity  An estimate of the atom density, used for peformance optimization and only with \p gridIndex = 0
472  * \param[in]     atomInfo     Atom information flags
473  * \param[in]     x            Coordinates for atoms to grid
474  * \param[in]     numAtomsMoved  The number of atoms that will move to another domain, pass 0 without DD
475  * \param[in]     move         Move flags for atoms, pass nullptr without DD
476  */
477 void nbnxn_put_on_grid(nonbonded_verlet_t             *nb_verlet,
478                        const matrix                    box,
479                        int                             gridIndex,
480                        const rvec                      lowerCorner,
481                        const rvec                      upperCorner,
482                        const gmx::UpdateGroupsCog     *updateGroupsCog,
483                        gmx::Range<int>                 atomRange,
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