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