bcb3b34c3379835c06a9529a3094f7f4f193e0e2
[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 Nbnxm
159 {
160 enum class KernelType;
161 }
162
163 namespace Nbnxm
164 {
165
166 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
167 enum class KernelType : int
168 {
169     NotSet = 0,
170     Cpu4x4_PlainC,
171     Cpu4xN_Simd_4xN,
172     Cpu4xN_Simd_2xNN,
173     Gpu8x8x8,
174     Cpu8x8x8_PlainC,
175     Count
176 };
177
178 /*! \brief Ewald exclusion types */
179 enum class EwaldExclusionType : int
180 {
181     NotSet = 0,
182     Table,
183     Analytical,
184     DecidedByGpuModule
185 };
186
187 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
188 struct KernelSetup
189 {
190     //! The non-bonded type, also affects the pairlist construction kernel
191     KernelType kernelType = KernelType::NotSet;
192     //! Ewald exclusion computation handling type, currently only used for CPU
193     EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
194 };
195
196 /*! \brief Return a string identifying the kernel type.
197  *
198  * \param [in] kernelType   nonbonded kernel type, takes values from the nbnxn_kernel_type enum
199  * \returns                 a string identifying the kernel corresponding to the type passed as argument
200  */
201 const char* lookup_kernel_name(Nbnxm::KernelType kernelType);
202
203 } // namespace Nbnxm
204
205 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
206 enum
207 {
208     enbvClearFNo,
209     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                        NbnxmGpu*                         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 { return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8; }
229
230     //! Returns whether a GPU is emulated for the non-bonded calculations
231     bool emulateGpu() const
232     {
233         return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
234     }
235
236     //! Return whether the pairlist is of simple, CPU type
237     bool pairlistIsSimple() const { return !useGpu() && !emulateGpu(); }
238
239     //! Initialize the pair list sets, TODO this should be private
240     void initPairlistSets(bool haveMultipleDomains);
241
242     //! Returns the order of the local atoms on the grid
243     gmx::ArrayRef<const int> getLocalAtomOrder() const;
244
245     //! Sets the order of the local atoms to the order grid atom ordering
246     void setLocalAtomOrder();
247
248     //! Returns the index position of the atoms on the search grid
249     gmx::ArrayRef<const int> getGridIndices() const;
250
251     /*! \brief Constructs the pairlist for the given locality
252      *
253      * When there are no non-self exclusions, \p exclusions can be empty.
254      * Otherwise the number of lists in \p exclusions should match the number
255      * of atoms when not using DD, or the total number of atoms in the i-zones
256      * when using DD.
257      *
258      * \param[in] iLocality   The interaction locality: local or non-local
259      * \param[in] exclusions  Lists of exclusions for every atom.
260      * \param[in] step        Used to set the list creation step
261      * \param[in,out] nrnb    Flop accounting struct, can be nullptr
262      */
263     void constructPairlist(gmx::InteractionLocality     iLocality,
264                            const gmx::ListOfLists<int>& exclusions,
265                            int64_t                      step,
266                            t_nrnb*                      nrnb);
267
268     //! Updates all the atom properties in Nbnxm
269     void setAtomProperties(gmx::ArrayRef<const int>  atomTypes,
270                            gmx::ArrayRef<const real> atomCharges,
271                            gmx::ArrayRef<const int>  atomInfo);
272
273     /*!\brief Convert the coordinates to NBNXM format for the given locality.
274      *
275      * The API function for the transformation of the coordinates from one layout to another.
276      *
277      * \param[in] locality     Whether coordinates for local or non-local atoms should be
278      * transformed. \param[in] fillLocal    If the coordinates for filler particles should be
279      * zeroed. \param[in] coordinates  Coordinates in plain rvec format to be transformed.
280      */
281     void convertCoordinates(gmx::AtomLocality locality, bool fillLocal, gmx::ArrayRef<const gmx::RVec> coordinates);
282
283     /*!\brief Convert the coordinates to NBNXM format on the GPU for the given locality
284      *
285      * The API function for the transformation of the coordinates from one layout to another in the GPU memory.
286      *
287      * \param[in] locality        Whether coordinates for local or non-local atoms should be transformed.
288      * \param[in] fillLocal       If the coordinates for filler particles should be zeroed.
289      * \param[in] d_x             GPU coordinates buffer in plain rvec format to be transformed.
290      * \param[in] xReadyOnDevice  Event synchronizer indicating that the coordinates are ready in the device memory.
291      */
292     void convertCoordinatesGpu(gmx::AtomLocality       locality,
293                                bool                    fillLocal,
294                                DeviceBuffer<gmx::RVec> d_x,
295                                GpuEventSynchronizer*   xReadyOnDevice);
296
297     //! Init for GPU version of setup coordinates in Nbnxm
298     void atomdata_init_copy_x_to_nbat_x_gpu();
299
300     //! Sync the nonlocal GPU stream with dependent tasks in the local queue.
301     void insertNonlocalGpuDependency(gmx::InteractionLocality interactionLocality);
302
303     //! Returns a reference to the pairlist sets
304     const PairlistSets& pairlistSets() const { return *pairlistSets_; }
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(gmx::InteractionLocality iLocality, const rvec* shift_vec);
314
315     //! Dispatches the dynamic pruning kernel for GPU lists
316     void dispatchPruneKernelGpu(int64_t step);
317
318     //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
319     void dispatchNonbondedKernel(gmx::InteractionLocality   iLocality,
320                                  const interaction_const_t& ic,
321                                  const gmx::StepWorkload&   stepWork,
322                                  int                        clearF,
323                                  const t_forcerec&          fr,
324                                  gmx_enerdata_t*            enerd,
325                                  t_nrnb*                    nrnb);
326
327     //! Executes the non-bonded free-energy kernel, always runs on the CPU
328     void dispatchFreeEnergyKernel(gmx::InteractionLocality   iLocality,
329                                   const t_forcerec*          fr,
330                                   rvec                       x[],
331                                   gmx::ForceWithShiftForces* forceWithShiftForces,
332                                   const t_mdatoms&           mdatoms,
333                                   t_lambda*                  fepvals,
334                                   gmx::ArrayRef<const real>  lambda,
335                                   gmx_enerdata_t*            enerd,
336                                   const gmx::StepWorkload&   stepWork,
337                                   t_nrnb*                    nrnb);
338
339     /*! \brief Add the forces stored in nbat to f, zeros the forces in nbat
340      * \param [in] locality         Local or non-local
341      * \param [inout] force         Force to be added to
342      */
343     void atomdata_add_nbat_f_to_f(gmx::AtomLocality locality, gmx::ArrayRef<gmx::RVec> force);
344
345     /*! \brief Add the forces stored in nbat to total force using GPU buffer opse
346      *
347      * \param [in]     locality             Local or non-local
348      * \param [in,out] totalForcesDevice    Force to be added to
349      * \param [in]     forcesPmeDevice      Device buffer with PME forces
350      * \param[in]      dependencyList       List of synchronizers that represent the dependencies the reduction task needs to sync on.
351      * \param [in]     useGpuFPmeReduction  Whether PME forces should be added
352      * \param [in]     accumulateForce      If the total force buffer already contains data
353      */
354     void atomdata_add_nbat_f_to_f_gpu(gmx::AtomLocality                          locality,
355                                       DeviceBuffer<gmx::RVec>                    totalForcesDevice,
356                                       void*                                      forcesPmeDevice,
357                                       gmx::ArrayRef<GpuEventSynchronizer* const> dependencyList,
358                                       bool useGpuFPmeReduction,
359                                       bool accumulateForce);
360
361     /*! \brief Outer body of function to perform initialization for F buffer operations on GPU.
362      *
363      * \param localReductionDone     Pointer to an event synchronizer that marks the completion of the local f buffer ops kernel.
364      */
365     void atomdata_init_add_nbat_f_to_f_gpu(GpuEventSynchronizer* localReductionDone);
366
367     /*! \brief return GPU pointer to f in rvec format */
368     void* get_gpu_frvec();
369
370     //! Return the kernel setup
371     const Nbnxm::KernelSetup& kernelSetup() const { return kernelSetup_; }
372
373     //! Returns the outer radius for the pair list
374     real pairlistInnerRadius() const;
375
376     //! Returns the outer radius for the pair list
377     real pairlistOuterRadius() const;
378
379     //! Changes the pair-list outer and inner radius
380     void changePairlistRadii(real rlistOuter, real rlistInner);
381
382     //! Set up internal flags that indicate what type of short-range work there is.
383     void setupGpuShortRangeWork(const gmx::GpuBonded* gpuBonded, gmx::InteractionLocality iLocality);
384
385     // TODO: Make all data members private
386 public:
387     //! All data related to the pair lists
388     std::unique_ptr<PairlistSets> pairlistSets_;
389     //! Working data for constructing the pairlists
390     std::unique_ptr<PairSearch> pairSearch_;
391     //! Atom data
392     std::unique_ptr<nbnxn_atomdata_t> nbat;
393
394 private:
395     //! The non-bonded setup, also affects the pairlist construction kernel
396     Nbnxm::KernelSetup kernelSetup_;
397     //! \brief Pointer to wallcycle structure.
398     gmx_wallcycle* wcycle_;
399
400 public:
401     //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
402     NbnxmGpu* gpu_nbv;
403 };
404
405 namespace Nbnxm
406 {
407
408 /*! \brief Creates an Nbnxm object */
409 std::unique_ptr<nonbonded_verlet_t> init_nb_verlet(const gmx::MDLogger& mdlog,
410                                                    const t_inputrec*    ir,
411                                                    const t_forcerec*    fr,
412                                                    const t_commrec*     cr,
413                                                    const gmx_hw_info_t& hardwareInfo,
414                                                    bool                 useGpuForNonbonded,
415                                                    const gmx::DeviceStreamManager* deviceStreamManager,
416                                                    const gmx_mtop_t*               mtop,
417                                                    matrix                          box,
418                                                    gmx_wallcycle*                  wcycle);
419
420 } // namespace Nbnxm
421
422 /*! \brief Put the atoms on the pair search grid.
423  *
424  * Only atoms with indices wihtin \p atomRange in x are put on the grid.
425  * When \p updateGroupsCog != nullptr, atoms are put on the grid
426  * based on the center of geometry of the group they belong to.
427  * Atoms or COGs of groups should be within the bounding box provided,
428  * this is checked in debug builds when not using update groups.
429  * The atom density is used to determine the grid size when \p gridIndex = 0.
430  * When \p atomDensity <= 0, the density is determined from atomEnd-atomStart
431  * and the bounding box corners.
432  * With domain decomposition, part of the atoms might have migrated,
433  * but have not been removed yet. This count is given by \p numAtomsMoved.
434  * When \p move[i] < 0 particle i has migrated and will not be put on the grid.
435  *
436  * \param[in,out] nb_verlet    The non-bonded object
437  * \param[in]     box          Box used for periodic distance calculations
438  * \param[in]     gridIndex    The index of the grid to spread to, always 0 except with test particle insertion
439  * \param[in]     lowerCorner  Atom groups to be gridded should have coordinates >= this corner
440  * \param[in]     upperCorner  Atom groups to be gridded should have coordinates <= this corner
441  * \param[in]     updateGroupsCog  Centers of geometry for update groups, pass nullptr when not using update groups
442  * \param[in]     atomRange    Range of atoms to grid
443  * \param[in]     atomDensity  An estimate of the atom density, used for peformance optimization and only with \p gridIndex = 0
444  * \param[in]     atomInfo     Atom information flags
445  * \param[in]     x            Coordinates for atoms to grid
446  * \param[in]     numAtomsMoved  The number of atoms that will move to another domain, pass 0 without DD
447  * \param[in]     move         Move flags for atoms, pass nullptr without DD
448  */
449 void nbnxn_put_on_grid(nonbonded_verlet_t*            nb_verlet,
450                        const matrix                   box,
451                        int                            gridIndex,
452                        const rvec                     lowerCorner,
453                        const rvec                     upperCorner,
454                        const gmx::UpdateGroupsCog*    updateGroupsCog,
455                        gmx::Range<int>                atomRange,
456                        real                           atomDensity,
457                        gmx::ArrayRef<const int>       atomInfo,
458                        gmx::ArrayRef<const gmx::RVec> x,
459                        int                            numAtomsMoved,
460                        const int*                     move);
461
462 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
463  *
464  * with domain decomposition. Should be called after calling
465  * nbnxn_search_put_on_grid for the local atoms / home zone.
466  */
467 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t*              nb_verlet,
468                                 const struct gmx_domdec_zones_t* zones,
469                                 gmx::ArrayRef<const int>         atomInfo,
470                                 gmx::ArrayRef<const gmx::RVec>   x);
471
472 #endif // GMX_NBNXN_NBNXM_H