Minor nbnxm cleanup
[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/nbnxm/pairlist.h"
105 #include "gromacs/utility/arrayref.h"
106 #include "gromacs/utility/enumerationhelpers.h"
107 #include "gromacs/utility/real.h"
108
109 #include "locality.h"
110
111 // TODO: Remove this include and the two nbnxm includes above
112 #include "nbnxm_gpu.h"
113
114 struct gmx_device_info_t;
115 struct gmx_domdec_zones_t;
116 struct gmx_enerdata_t;
117 struct gmx_hw_info_t;
118 struct gmx_mtop_t;
119 struct interaction_const_t;
120 struct nbnxn_pairlist_set_t;
121 struct nbnxn_search;
122 struct nonbonded_verlet_t;
123 struct t_blocka;
124 struct t_commrec;
125 struct t_nrnb;
126 struct t_forcerec;
127 struct t_inputrec;
128
129 namespace gmx
130 {
131 class MDLogger;
132 class UpdateGroupsCog;
133 }
134
135 /*! \brief Resources that can be used to execute non-bonded kernels on */
136 enum class NonbondedResource : int
137 {
138     Cpu,
139     Gpu,
140     EmulateGpu
141 };
142
143 namespace Nbnxm
144 {
145
146 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
147 enum class KernelType : int
148 {
149     NotSet = 0,
150     Cpu4x4_PlainC,
151     Cpu4xN_Simd_4xN,
152     Cpu4xN_Simd_2xNN,
153     Gpu8x8x8,
154     Cpu8x8x8_PlainC,
155     Count
156 };
157
158 /*! \brief Ewald exclusion types */
159 enum class EwaldExclusionType : int
160 {
161     NotSet = 0,
162     Table,
163     Analytical,
164     DecidedByGpuModule
165 };
166
167 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
168 struct KernelSetup
169 {
170     //! The non-bonded type, also affects the pairlist construction kernel
171     KernelType         kernelType = KernelType::NotSet;
172     //! Ewald exclusion computation handling type, currently only used for CPU
173     EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
174 };
175
176 /*! \brief Return a string identifying the kernel type.
177  *
178  * \param [in] kernelType   nonbonded kernel type, takes values from the nbnxn_kernel_type enum
179  * \returns                 a string identifying the kernel corresponding to the type passed as argument
180  */
181 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
182
183 } // namespace Nbnxm
184
185 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
186 enum {
187     enbvClearFNo, enbvClearFYes
188 };
189
190 /*! \brief Generates a pair-list for the given locality.
191  *
192  * With perturbed particles, also a group scheme style nbl_fep list is made.
193  */
194 void nbnxn_make_pairlist(nonbonded_verlet_t         *nbv,
195                          Nbnxm::InteractionLocality  iLocality,
196                          nbnxn_pairlist_set_t       *pairlistSet,
197                          const t_blocka             *excl,
198                          int64_t                     step,
199                          t_nrnb                     *nrnb);
200
201 /*! \brief Prune all pair-lists with given locality (currently CPU only)
202  *
203  * For all pair-lists with given locality, takes the outer list and prunes out
204  * pairs beyond the pairlist inner radius and writes the result to a list that is
205  * to be consumed by the non-bonded kernel.
206  */
207 void NbnxnDispatchPruneKernel(nbnxn_pairlist_set_t   *pairlistSet,
208                               Nbnxm::KernelType       kernelType,
209                               const nbnxn_atomdata_t *nbat,
210                               const rvec             *shift_vec);
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         //! Returns whether a GPU is use for the non-bonded calculations
218         bool useGpu() const
219         {
220             return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
221         }
222
223         //! Returns whether a GPU is emulated for the non-bonded calculations
224         bool emulateGpu() const
225         {
226             return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
227         }
228
229         //! Return whether the pairlist is of simple, CPU type
230         bool pairlistIsSimple() const
231         {
232             return !useGpu() && !emulateGpu();
233         }
234
235         //! Initialize the pair list sets, TODO this should be private
236         void initPairlistSets(bool haveMultipleDomains);
237
238         //! Returns a reference to the pairlist set for the requested locality
239         const nbnxn_pairlist_set_t &pairlistSet(Nbnxm::InteractionLocality iLocality) const
240         {
241             GMX_ASSERT(static_cast<size_t>(iLocality) < pairlistSets_.size(),
242                        "The requested locality should be in the list");
243             return pairlistSets_[static_cast<int>(iLocality)];
244         }
245
246         //! Constructs the pairlist for the given locality
247         void constructPairlist(Nbnxm::InteractionLocality  iLocality,
248                                const t_blocka             *excl,
249                                int64_t                     step,
250                                t_nrnb                     *nrnb)
251         {
252             nbnxn_make_pairlist(this, iLocality, &pairlistSets_[static_cast<int>(iLocality)], excl, step, nrnb);
253         }
254
255         //! Dispatches the dynamic pruning kernel for the given locality
256         void dispatchPruneKernel(Nbnxm::InteractionLocality  iLocality,
257                                  const rvec                 *shift_vec)
258         {
259             GMX_ASSERT(static_cast<size_t>(iLocality) < pairlistSets_.size(),
260                        "The requested locality should be in the list");
261             NbnxnDispatchPruneKernel(&pairlistSets_[static_cast<int>(iLocality)],
262                                      kernelSetup_.kernelType, nbat, shift_vec);
263         }
264
265         //! Return the kernel setup
266         const Nbnxm::KernelSetup &kernelSetup() const
267         {
268             return kernelSetup_;
269         }
270
271         //! Sets the kernel setup, TODO: make private
272         void setKernelSetup(const Nbnxm::KernelSetup &kernelSetup)
273         {
274             kernelSetup_ = kernelSetup;
275         }
276
277         //! Returns the a list of free-energy pairlists for the given locality
278         const gmx::ArrayRef<t_nblist const * const>
279         freeEnergyPairlistSet(Nbnxm::InteractionLocality iLocality) const
280         {
281             return pairlistSet(iLocality).nbl_fep;
282         }
283
284         //! Parameters for the search and list pruning setup
285         std::unique_ptr<NbnxnListParameters>  listParams;
286         //! Working data for constructing the pairlists
287         std::unique_ptr<nbnxn_search>         nbs;
288     private:
289         //! Local and, optionally, non-local pairlist sets
290         std::vector<nbnxn_pairlist_set_t>     pairlistSets_;
291     public:
292         //! Atom data
293         nbnxn_atomdata_t                     *nbat;
294
295     private:
296         //! The non-bonded setup, also affects the pairlist construction kernel
297         Nbnxm::KernelSetup   kernelSetup_;
298     public:
299
300         gmx_nbnxn_gpu_t     *gpu_nbv;         /**< pointer to GPU nb verlet data     */
301         int                  min_ci_balanced; /**< pair list balancing parameter used for the 8x8x8 GPU kernels    */
302 };
303
304 namespace Nbnxm
305 {
306
307 /*! \brief Initializes the nbnxn module */
308 void init_nb_verlet(const gmx::MDLogger     &mdlog,
309                     nonbonded_verlet_t     **nb_verlet,
310                     gmx_bool                 bFEP_NonBonded,
311                     const t_inputrec        *ir,
312                     const t_forcerec        *fr,
313                     const t_commrec         *cr,
314                     const gmx_hw_info_t     &hardwareInfo,
315                     const gmx_device_info_t *deviceInfo,
316                     const gmx_mtop_t        *mtop,
317                     matrix                   box);
318
319 } // namespace Nbnxm
320
321 /*! \brief Put the atoms on the pair search grid.
322  *
323  * Only atoms atomStart to atomEnd in x are put on the grid.
324  * The atom_density is used to determine the grid size.
325  * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
326  * With domain decomposition part of the n particles might have migrated,
327  * but have not been removed yet. This count is given by nmoved.
328  * When move[i] < 0 particle i has migrated and will not be put on the grid.
329  * Without domain decomposition move will be NULL.
330  */
331 void nbnxn_put_on_grid(nonbonded_verlet_t             *nb_verlet,
332                        const matrix                    box,
333                        int                             ddZone,
334                        const rvec                      lowerCorner,
335                        const rvec                      upperCorner,
336                        const gmx::UpdateGroupsCog     *updateGroupsCog,
337                        int                             atomStart,
338                        int                             atomEnd,
339                        real                            atomDensity,
340                        const int                      *atinfo,
341                        gmx::ArrayRef<const gmx::RVec>  x,
342                        int                             numAtomsMoved,
343                        const int                      *move);
344
345 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
346  *
347  * with domain decomposition. Should be called after calling
348  * nbnxn_search_put_on_grid for the local atoms / home zone.
349  */
350 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t              *nb_verlet,
351                                 const struct gmx_domdec_zones_t *zones,
352                                 const int                       *atinfo,
353                                 gmx::ArrayRef<const gmx::RVec>   x);
354
355 /*! \brief Returns the number of x and y cells in the local grid */
356 void nbnxn_get_ncells(const nbnxn_search *nbs, int *ncx, int *ncy);
357
358 /*! \brief Returns the order indices of the atoms on the pairlist search grid */
359 gmx::ArrayRef<const int> nbnxn_get_atomorder(const nbnxn_search* nbs);
360
361 /*! \brief Renumbers the atom indices on the grid to consecutive order */
362 void nbnxn_set_atomorder(nbnxn_search *nbs);
363
364 /*! \brief Returns the index position of the atoms on the pairlist search grid */
365 gmx::ArrayRef<const int> nbnxn_get_gridindices(const nbnxn_search* nbs);
366
367 /*! \brief Returns the number of steps performed with the current pair list */
368 int nbnxnNumStepsWithPairlist(const nonbonded_verlet_t   &nbv,
369                               Nbnxm::InteractionLocality  ilocality,
370                               int64_t                     step);
371
372 /*! \brief Returns whether step is a dynamic list pruning step */
373 bool nbnxnIsDynamicPairlistPruningStep(const nonbonded_verlet_t   &nbv,
374                                        Nbnxm::InteractionLocality  ilocality,
375                                        int64_t                     step);
376
377 /*! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU */
378 void NbnxnDispatchKernel(nonbonded_verlet_t         *nbv,
379                          Nbnxm::InteractionLocality  iLocality,
380                          const interaction_const_t  &ic,
381                          int                         forceFlags,
382                          int                         clearF,
383                          t_forcerec                 *fr,
384                          gmx_enerdata_t             *enerd,
385                          t_nrnb                     *nrnb);
386
387 #endif // GMX_NBNXN_NBNXN_H