Apply clang-format to source tree
[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/mdtypes/locality.h"
118 #include "gromacs/utility/arrayref.h"
119 #include "gromacs/utility/enumerationhelpers.h"
120 #include "gromacs/utility/range.h"
121 #include "gromacs/utility/real.h"
122
123 // TODO: Remove this include
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 } // namespace gmx
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 {
210     enbvClearFNo,
211     enbvClearFYes
212 };
213
214 /*! \libinternal
215  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
216 struct nonbonded_verlet_t
217 {
218 public:
219     //! Constructs an object from its components
220     nonbonded_verlet_t(std::unique_ptr<PairlistSets>     pairlistSets,
221                        std::unique_ptr<PairSearch>       pairSearch,
222                        std::unique_ptr<nbnxn_atomdata_t> nbat,
223                        const Nbnxm::KernelSetup&         kernelSetup,
224                        gmx_nbnxn_gpu_t*                  gpu_nbv,
225                        gmx_wallcycle*                    wcycle);
226
227     ~nonbonded_verlet_t();
228
229     //! Returns whether a GPU is use for the non-bonded calculations
230     bool useGpu() const { return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8; }
231
232     //! Returns whether a GPU is emulated for the non-bonded calculations
233     bool emulateGpu() const
234     {
235         return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
236     }
237
238     //! Return whether the pairlist is of simple, CPU type
239     bool pairlistIsSimple() const { return !useGpu() && !emulateGpu(); }
240
241     //! Initialize the pair list sets, TODO this should be private
242     void initPairlistSets(bool haveMultipleDomains);
243
244     //! Returns the order of the local atoms on the grid
245     gmx::ArrayRef<const int> getLocalAtomOrder() const;
246
247     //! Sets the order of the local atoms to the order grid atom ordering
248     void setLocalAtomOrder();
249
250     //! Returns the index position of the atoms on the search grid
251     gmx::ArrayRef<const int> getGridIndices() const;
252
253     //! Constructs the pairlist for the given locality
254     void constructPairlist(gmx::InteractionLocality iLocality, const t_blocka* excl, int64_t step, t_nrnb* nrnb);
255
256     //! Updates all the atom properties in Nbnxm
257     void setAtomProperties(const t_mdatoms& mdatoms, gmx::ArrayRef<const int> atomInfo);
258
259     /*!\brief Convert the coordinates to NBNXM format for the given locality.
260      *
261      * The API function for the transformation of the coordinates from one layout to another.
262      *
263      * \param[in] locality     Whether coordinates for local or non-local atoms should be
264      * transformed. \param[in] fillLocal    If the coordinates for filler particles should be
265      * zeroed. \param[in] coordinates  Coordinates in plain rvec format to be transformed.
266      */
267     void convertCoordinates(gmx::AtomLocality locality, bool fillLocal, gmx::ArrayRef<const gmx::RVec> coordinates);
268
269     /*!\brief Convert the coordinates to NBNXM format on the GPU for the given locality
270      *
271      * The API function for the transformation of the coordinates from one layout to another in the GPU memory.
272      *
273      * \param[in] locality        Whether coordinates for local or non-local atoms should be transformed.
274      * \param[in] fillLocal       If the coordinates for filler particles should be zeroed.
275      * \param[in] d_x             GPU coordinates buffer in plain rvec format to be transformed.
276      * \param[in] xReadyOnDevice  Event synchronizer indicating that the coordinates are ready in the device memory.
277      */
278     void convertCoordinatesGpu(gmx::AtomLocality     locality,
279                                bool                  fillLocal,
280                                DeviceBuffer<float>   d_x,
281                                GpuEventSynchronizer* xReadyOnDevice);
282
283     //! Init for GPU version of setup coordinates in Nbnxm
284     void atomdata_init_copy_x_to_nbat_x_gpu();
285
286     //! Sync the nonlocal GPU stream with dependent tasks in the local queue.
287     void insertNonlocalGpuDependency(gmx::InteractionLocality interactionLocality);
288
289     //! Returns a reference to the pairlist sets
290     const PairlistSets& pairlistSets() const { return *pairlistSets_; }
291
292     //! Returns whether step is a dynamic list pruning step, for CPU lists
293     bool isDynamicPruningStepCpu(int64_t step) const;
294
295     //! Returns whether step is a dynamic list pruning step, for GPU lists
296     bool isDynamicPruningStepGpu(int64_t step) const;
297
298     //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
299     void dispatchPruneKernelCpu(gmx::InteractionLocality iLocality, const rvec* shift_vec);
300
301     //! Dispatches the dynamic pruning kernel for GPU lists
302     void dispatchPruneKernelGpu(int64_t step);
303
304     //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
305     void dispatchNonbondedKernel(gmx::InteractionLocality   iLocality,
306                                  const interaction_const_t& ic,
307                                  const gmx::StepWorkload&   stepWork,
308                                  int                        clearF,
309                                  const t_forcerec&          fr,
310                                  gmx_enerdata_t*            enerd,
311                                  t_nrnb*                    nrnb);
312
313     //! Executes the non-bonded free-energy kernel, always runs on the CPU
314     void dispatchFreeEnergyKernel(gmx::InteractionLocality   iLocality,
315                                   const t_forcerec*          fr,
316                                   rvec                       x[],
317                                   gmx::ForceWithShiftForces* forceWithShiftForces,
318                                   const t_mdatoms&           mdatoms,
319                                   t_lambda*                  fepvals,
320                                   real*                      lambda,
321                                   gmx_enerdata_t*            enerd,
322                                   const gmx::StepWorkload&   stepWork,
323                                   t_nrnb*                    nrnb);
324
325     /*! \brief Add the forces stored in nbat to f, zeros the forces in nbat
326      * \param [in] locality         Local or non-local
327      * \param [inout] force         Force to be added to
328      */
329     void atomdata_add_nbat_f_to_f(gmx::AtomLocality locality, gmx::ArrayRef<gmx::RVec> force);
330
331     /*! \brief Add the forces stored in nbat to total force using GPU buffer opse
332      *
333      * \param [in]     locality             Local or non-local
334      * \param [in,out] totalForcesDevice    Force to be added to
335      * \param [in]     forcesPmeDevice      Device buffer with PME forces
336      * \param[in]      dependencyList       List of synchronizers that represent the dependencies the reduction task needs to sync on.
337      * \param [in]     useGpuFPmeReduction  Whether PME forces should be added
338      * \param [in]     accumulateForce      If the total force buffer already contains data
339      */
340     void atomdata_add_nbat_f_to_f_gpu(gmx::AtomLocality                          locality,
341                                       DeviceBuffer<float>                        totalForcesDevice,
342                                       void*                                      forcesPmeDevice,
343                                       gmx::ArrayRef<GpuEventSynchronizer* const> dependencyList,
344                                       bool useGpuFPmeReduction,
345                                       bool accumulateForce);
346
347     /*! \brief Outer body of function to perform initialization for F buffer operations on GPU.
348      *
349      * \param localReductionDone     Pointer to an event synchronizer that marks the completion of the local f buffer ops kernel.
350      */
351     void atomdata_init_add_nbat_f_to_f_gpu(GpuEventSynchronizer* localReductionDone);
352
353     /*! \brief Wait for non-local copy of coordinate buffer from device to host */
354     void wait_nonlocal_x_copy_D2H_done();
355
356     /*! \brief return GPU pointer to f in rvec format */
357     void* get_gpu_frvec();
358
359     /*! \brief Ensure local stream waits for non-local stream */
360     void stream_local_wait_for_nonlocal();
361
362     //! Return the kernel setup
363     const Nbnxm::KernelSetup& kernelSetup() const { return kernelSetup_; }
364
365     //! Returns the outer radius for the pair list
366     real pairlistInnerRadius() const;
367
368     //! Returns the outer radius for the pair list
369     real pairlistOuterRadius() const;
370
371     //! Changes the pair-list outer and inner radius
372     void changePairlistRadii(real rlistOuter, real rlistInner);
373
374     //! Set up internal flags that indicate what type of short-range work there is.
375     void setupGpuShortRangeWork(const gmx::GpuBonded* gpuBonded, const gmx::InteractionLocality iLocality)
376     {
377         if (useGpu() && !emulateGpu())
378         {
379             Nbnxm::setupGpuShortRangeWork(gpu_nbv, gpuBonded, iLocality);
380         }
381     }
382
383     //! Returns true if there is GPU short-range work for the given atom locality.
384     bool haveGpuShortRangeWork(const gmx::AtomLocality aLocality)
385     {
386         return ((useGpu() && !emulateGpu()) && Nbnxm::haveGpuShortRangeWork(gpu_nbv, aLocality));
387     }
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_nbnxn_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