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