Separate CPU NB kernel and buffer clearing subcounters
[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 \file
84  *
85  * \brief This file contains the public interface of the nbnxm module
86  * that implements the NxM atom cluster non-bonded algorithm to efficiently
87  * compute pair forces.
88  *
89  *
90  * \author Berk Hess <hess@kth.se>
91  * \author Szilárd Páll <pall.szilard@gmail.com>
92  *
93  * \inlibraryapi
94  * \ingroup __module_nbnxm
95  */
96
97
98 #ifndef GMX_NBNXM_NBNXM_H
99 #define GMX_NBNXM_NBNXM_H
100
101 #include <memory>
102
103 #include "gromacs/math/vectypes.h"
104 #include "gromacs/utility/arrayref.h"
105 #include "gromacs/utility/enumerationhelpers.h"
106 #include "gromacs/utility/real.h"
107
108 #include "locality.h"
109
110 // TODO: Remove this include and the two nbnxm includes above
111 #include "nbnxm_gpu.h"
112
113 struct gmx_device_info_t;
114 struct gmx_domdec_zones_t;
115 struct gmx_enerdata_t;
116 struct gmx_hw_info_t;
117 struct gmx_mtop_t;
118 struct gmx_wallcycle;
119 struct interaction_const_t;
120 struct nonbonded_verlet_t;
121 class PairSearch;
122 class PairlistSets;
123 struct t_blocka;
124 struct t_commrec;
125 struct t_lambda;
126 struct t_mdatoms;
127 struct t_nrnb;
128 struct t_forcerec;
129 struct t_inputrec;
130
131 namespace gmx
132 {
133 class MDLogger;
134 class UpdateGroupsCog;
135 }
136
137 namespace Nbnxm
138 {
139 enum class KernelType;
140 }
141
142 namespace Nbnxm
143 {
144
145 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
146 enum class KernelType : int
147 {
148     NotSet = 0,
149     Cpu4x4_PlainC,
150     Cpu4xN_Simd_4xN,
151     Cpu4xN_Simd_2xNN,
152     Gpu8x8x8,
153     Cpu8x8x8_PlainC,
154     Count
155 };
156
157 /*! \brief Ewald exclusion types */
158 enum class EwaldExclusionType : int
159 {
160     NotSet = 0,
161     Table,
162     Analytical,
163     DecidedByGpuModule
164 };
165
166 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
167 struct KernelSetup
168 {
169     //! The non-bonded type, also affects the pairlist construction kernel
170     KernelType         kernelType = KernelType::NotSet;
171     //! Ewald exclusion computation handling type, currently only used for CPU
172     EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
173 };
174
175 /*! \brief Return a string identifying the kernel type.
176  *
177  * \param [in] kernelType   nonbonded kernel type, takes values from the nbnxn_kernel_type enum
178  * \returns                 a string identifying the kernel corresponding to the type passed as argument
179  */
180 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
181
182 } // namespace Nbnxm
183
184 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
185 enum {
186     enbvClearFNo, enbvClearFYes
187 };
188
189 /*! \libinternal
190  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
191 struct nonbonded_verlet_t
192 {
193     public:
194         //! Constructs an object from its components
195         nonbonded_verlet_t(std::unique_ptr<PairlistSets>      pairlistSets,
196                            std::unique_ptr<PairSearch>        pairSearch,
197                            std::unique_ptr<nbnxn_atomdata_t>  nbat,
198                            const Nbnxm::KernelSetup          &kernelSetup,
199                            gmx_nbnxn_gpu_t                   *gpu_nbv);
200
201         ~nonbonded_verlet_t();
202
203         //! Returns whether a GPU is use for the non-bonded calculations
204         bool useGpu() const
205         {
206             return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
207         }
208
209         //! Returns whether a GPU is emulated for the non-bonded calculations
210         bool emulateGpu() const
211         {
212             return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
213         }
214
215         //! Return whether the pairlist is of simple, CPU type
216         bool pairlistIsSimple() const
217         {
218             return !useGpu() && !emulateGpu();
219         }
220
221         //! Initialize the pair list sets, TODO this should be private
222         void initPairlistSets(bool haveMultipleDomains);
223
224         //! Returns the order of the local atoms on the grid
225         gmx::ArrayRef<const int> getLocalAtomOrder() const;
226
227         //! Sets the order of the local atoms to the order grid atom ordering
228         void setLocalAtomOrder();
229
230         //! Returns the index position of the atoms on the search grid
231         gmx::ArrayRef<const int> getGridIndices() const;
232
233         //! Constructs the pairlist for the given locality
234         void constructPairlist(Nbnxm::InteractionLocality  iLocality,
235                                const t_blocka             *excl,
236                                int64_t                     step,
237                                t_nrnb                     *nrnb);
238
239         //! Updates all the atom properties in Nbnxm
240         void setAtomProperties(const t_mdatoms          &mdatoms,
241                                gmx::ArrayRef<const int>  atomInfo);
242
243         //! Updates the coordinates in Nbnxm for the given locality
244         void setCoordinates(Nbnxm::AtomLocality             locality,
245                             bool                            fillLocal,
246                             gmx::ArrayRef<const gmx::RVec>  x,
247                             bool                            useGpu,
248                             void                           *xPmeDevicePtr,
249                             gmx_wallcycle                  *wcycle);
250
251         //! Init for GPU version of setup coordinates in Nbnxm, for the given locality
252         void atomdata_init_copy_x_to_nbat_x_gpu(Nbnxm::AtomLocality        locality);
253
254
255         //! Returns a reference to the pairlist sets
256         const PairlistSets &pairlistSets() const
257         {
258             return *pairlistSets_;
259         }
260
261         //! Returns whether step is a dynamic list pruning step, for CPU lists
262         bool isDynamicPruningStepCpu(int64_t step) const;
263
264         //! Returns whether step is a dynamic list pruning step, for GPU lists
265         bool isDynamicPruningStepGpu(int64_t step) const;
266
267         //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
268         void dispatchPruneKernelCpu(Nbnxm::InteractionLocality  iLocality,
269                                     const rvec                 *shift_vec);
270
271         //! Dispatches the dynamic pruning kernel for GPU lists
272         void dispatchPruneKernelGpu(int64_t step);
273
274         //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
275         void dispatchNonbondedKernel(Nbnxm::InteractionLocality  iLocality,
276                                      const interaction_const_t  &ic,
277                                      int                         forceFlags,
278                                      int                         clearF,
279                                      const t_forcerec           &fr,
280                                      gmx_enerdata_t             *enerd,
281                                      t_nrnb                     *nrnb,
282                                      gmx_wallcycle              *wcycle);
283
284         //! Executes the non-bonded free-energy kernel, always runs on the CPU
285         void dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocality,
286                                       t_forcerec                 *fr,
287                                       rvec                        x[],
288                                       rvec                        f[],
289                                       const t_mdatoms            &mdatoms,
290                                       t_lambda                   *fepvals,
291                                       real                       *lambda,
292                                       gmx_enerdata_t             *enerd,
293                                       int                         forceFlags,
294                                       t_nrnb                     *nrnb);
295
296         //! Add the forces stored in nbat to f, zeros the forces in nbat */
297         void atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality  locality,
298                                       rvec                *f,
299                                       gmx_wallcycle       *wcycle);
300
301         //! Return the kernel setup
302         const Nbnxm::KernelSetup &kernelSetup() const
303         {
304             return kernelSetup_;
305         }
306
307         //! Returns the outer radius for the pair list
308         real pairlistInnerRadius() const;
309
310         //! Returns the outer radius for the pair list
311         real pairlistOuterRadius() const;
312
313         //! Changes the pair-list outer and inner radius
314         void changePairlistRadii(real rlistOuter,
315                                  real rlistInner);
316
317         // TODO: Make all data members private
318     public:
319         //! All data related to the pair lists
320         std::unique_ptr<PairlistSets>     pairlistSets_;
321         //! Working data for constructing the pairlists
322         std::unique_ptr<PairSearch>       pairSearch_;
323         //! Atom data
324         std::unique_ptr<nbnxn_atomdata_t> nbat;
325     private:
326         //! The non-bonded setup, also affects the pairlist construction kernel
327         Nbnxm::KernelSetup                kernelSetup_;
328     public:
329         //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
330         gmx_nbnxn_gpu_t                  *gpu_nbv;
331 };
332
333 namespace Nbnxm
334 {
335
336 /*! \brief Creates an Nbnxm object */
337 std::unique_ptr<nonbonded_verlet_t>
338 init_nb_verlet(const gmx::MDLogger     &mdlog,
339                gmx_bool                 bFEP_NonBonded,
340                const t_inputrec        *ir,
341                const t_forcerec        *fr,
342                const t_commrec         *cr,
343                const gmx_hw_info_t     &hardwareInfo,
344                const gmx_device_info_t *deviceInfo,
345                const gmx_mtop_t        *mtop,
346                matrix                   box);
347
348 } // namespace Nbnxm
349
350 /*! \brief Put the atoms on the pair search grid.
351  *
352  * Only atoms atomStart to atomEnd in x are put on the grid.
353  * The atom_density is used to determine the grid size.
354  * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
355  * With domain decomposition part of the n particles might have migrated,
356  * but have not been removed yet. This count is given by nmoved.
357  * When move[i] < 0 particle i has migrated and will not be put on the grid.
358  * Without domain decomposition move will be NULL.
359  */
360 void nbnxn_put_on_grid(nonbonded_verlet_t             *nb_verlet,
361                        const matrix                    box,
362                        int                             ddZone,
363                        const rvec                      lowerCorner,
364                        const rvec                      upperCorner,
365                        const gmx::UpdateGroupsCog     *updateGroupsCog,
366                        int                             atomStart,
367                        int                             atomEnd,
368                        real                            atomDensity,
369                        gmx::ArrayRef<const int>        atomInfo,
370                        gmx::ArrayRef<const gmx::RVec>  x,
371                        int                             numAtomsMoved,
372                        const int                      *move);
373
374 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
375  *
376  * with domain decomposition. Should be called after calling
377  * nbnxn_search_put_on_grid for the local atoms / home zone.
378  */
379 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t              *nb_verlet,
380                                 const struct gmx_domdec_zones_t *zones,
381                                 gmx::ArrayRef<const int>         atomInfo,
382                                 gmx::ArrayRef<const gmx::RVec>   x);
383
384 #endif // GMX_NBNXN_NBNXM_H