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