2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
36 // FIXME: remove the "__" prefix in front of the group def when we move the
37 // nonbonded code into separate dir.
39 /*! \libinternal \defgroup __module_nbnxm Short-range non-bonded interaction module
40 * \ingroup group_mdrun
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.
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.
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.
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.
70 * TODO: add a summary list and brief descriptions of the different submodules:
71 * search, CPU kernels, GPU glue code + kernels.
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>
80 * TODO: add more authors!
83 /*! \libinternal \file
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.
90 * \author Berk Hess <hess@kth.se>
91 * \author Szilárd Páll <pall.szilard@gmail.com>
94 * \ingroup __module_nbnxm
98 #ifndef GMX_NBNXM_NBNXM_H
99 #define GMX_NBNXM_NBNXM_H
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"
109 #include "locality.h"
111 // TODO: Remove this include and the two nbnxm includes above
112 #include "nbnxm_gpu.h"
114 struct gmx_device_info_t;
115 struct gmx_domdec_zones_t;
116 struct gmx_enerdata_t;
117 struct gmx_hw_info_t;
119 struct interaction_const_t;
120 struct nbnxn_pairlist_set_t;
121 struct nonbonded_verlet_t;
131 class UpdateGroupsCog;
134 /*! \brief Resources that can be used to execute non-bonded kernels on */
135 enum class NonbondedResource : int
145 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
146 enum class KernelType : int
157 /*! \brief Ewald exclusion types */
158 enum class EwaldExclusionType : int
166 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
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;
175 /*! \brief Return a string identifying the kernel type.
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
180 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
184 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
186 enbvClearFNo, enbvClearFYes
189 /*! \brief Generates a pair-list for the given locality.
191 * With perturbed particles, also a group scheme style nbl_fep list is made.
193 void nbnxn_make_pairlist(nonbonded_verlet_t *nbv,
194 Nbnxm::InteractionLocality iLocality,
195 nbnxn_pairlist_set_t *pairlistSet,
196 const t_blocka *excl,
200 /*! \brief Prune all pair-lists with given locality (currently CPU only)
202 * For all pair-lists with given locality, takes the outer list and prunes out
203 * pairs beyond the pairlist inner radius and writes the result to a list that is
204 * to be consumed by the non-bonded kernel.
206 void NbnxnDispatchPruneKernel(nbnxn_pairlist_set_t *pairlistSet,
207 Nbnxm::KernelType kernelType,
208 const nbnxn_atomdata_t *nbat,
209 const rvec *shift_vec);
212 * \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
213 struct nonbonded_verlet_t
216 //! Returns whether a GPU is use for the non-bonded calculations
219 return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
222 //! Returns whether a GPU is emulated for the non-bonded calculations
223 bool emulateGpu() const
225 return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
228 //! Return whether the pairlist is of simple, CPU type
229 bool pairlistIsSimple() const
231 return !useGpu() && !emulateGpu();
234 //! Initialize the pair list sets, TODO this should be private
235 void initPairlistSets(bool haveMultipleDomains);
237 //! Returns a reference to the pairlist set for the requested locality
238 const nbnxn_pairlist_set_t &pairlistSet(Nbnxm::InteractionLocality iLocality) const
240 GMX_ASSERT(static_cast<size_t>(iLocality) < pairlistSets_.size(),
241 "The requested locality should be in the list");
242 return pairlistSets_[static_cast<int>(iLocality)];
245 //! Constructs the pairlist for the given locality
246 void constructPairlist(Nbnxm::InteractionLocality iLocality,
247 const t_blocka *excl,
251 nbnxn_make_pairlist(this, iLocality, &pairlistSets_[static_cast<int>(iLocality)], excl, step, nrnb);
254 //! Dispatches the dynamic pruning kernel for the given locality
255 void dispatchPruneKernel(Nbnxm::InteractionLocality iLocality,
256 const rvec *shift_vec)
258 GMX_ASSERT(static_cast<size_t>(iLocality) < pairlistSets_.size(),
259 "The requested locality should be in the list");
260 NbnxnDispatchPruneKernel(&pairlistSets_[static_cast<int>(iLocality)],
261 kernelSetup_.kernelType, nbat, shift_vec);
264 //! Return the kernel setup
265 const Nbnxm::KernelSetup &kernelSetup() const
270 //! Sets the kernel setup, TODO: make private
271 void setKernelSetup(const Nbnxm::KernelSetup &kernelSetup)
273 kernelSetup_ = kernelSetup;
276 //! Returns the a list of free-energy pairlists for the given locality
277 const gmx::ArrayRef<t_nblist const * const>
278 freeEnergyPairlistSet(Nbnxm::InteractionLocality iLocality) const
280 return pairlistSet(iLocality).nbl_fep;
283 //! Parameters for the search and list pruning setup
284 std::unique_ptr<NbnxnListParameters> listParams;
285 //! Working data for constructing the pairlists
286 std::unique_ptr<nbnxn_search> nbs;
288 //! Local and, optionally, non-local pairlist sets
289 std::vector<nbnxn_pairlist_set_t> pairlistSets_;
292 nbnxn_atomdata_t *nbat;
295 //! The non-bonded setup, also affects the pairlist construction kernel
296 Nbnxm::KernelSetup kernelSetup_;
299 gmx_nbnxn_gpu_t *gpu_nbv; /**< pointer to GPU nb verlet data */
300 int min_ci_balanced; /**< pair list balancing parameter used for the 8x8x8 GPU kernels */
306 /*! \brief Initializes the nbnxn module */
307 void init_nb_verlet(const gmx::MDLogger &mdlog,
308 nonbonded_verlet_t **nb_verlet,
309 gmx_bool bFEP_NonBonded,
310 const t_inputrec *ir,
311 const t_forcerec *fr,
313 const gmx_hw_info_t &hardwareInfo,
314 const gmx_device_info_t *deviceInfo,
315 const gmx_mtop_t *mtop,
320 /*! \brief Put the atoms on the pair search grid.
322 * Only atoms atomStart to atomEnd in x are put on the grid.
323 * The atom_density is used to determine the grid size.
324 * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
325 * With domain decomposition part of the n particles might have migrated,
326 * but have not been removed yet. This count is given by nmoved.
327 * When move[i] < 0 particle i has migrated and will not be put on the grid.
328 * Without domain decomposition move will be NULL.
330 void nbnxn_put_on_grid(nonbonded_verlet_t *nb_verlet,
333 const rvec lowerCorner,
334 const rvec upperCorner,
335 const gmx::UpdateGroupsCog *updateGroupsCog,
340 gmx::ArrayRef<const gmx::RVec> x,
344 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
346 * with domain decomposition. Should be called after calling
347 * nbnxn_search_put_on_grid for the local atoms / home zone.
349 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t *nb_verlet,
350 const struct gmx_domdec_zones_t *zones,
352 gmx::ArrayRef<const gmx::RVec> x);
354 /*! \brief Returns the number of x and y cells in the local grid */
355 void nbnxn_get_ncells(nbnxn_search_t nbs, int *ncx, int *ncy);
357 /*! \brief Returns the order indices of the atoms on the pairlist search grid */
358 gmx::ArrayRef<const int> nbnxn_get_atomorder(const nbnxn_search* nbs);
360 /*! \brief Renumbers the atom indices on the grid to consecutive order */
361 void nbnxn_set_atomorder(nbnxn_search_t nbs);
363 /*! \brief Returns the index position of the atoms on the pairlist search grid */
364 gmx::ArrayRef<const int> nbnxn_get_gridindices(const nbnxn_search* nbs);
366 /*! \brief Returns the number of steps performed with the current pair list */
367 int nbnxnNumStepsWithPairlist(const nonbonded_verlet_t &nbv,
368 Nbnxm::InteractionLocality ilocality,
371 /*! \brief Returns whether step is a dynamic list pruning step */
372 bool nbnxnIsDynamicPairlistPruningStep(const nonbonded_verlet_t &nbv,
373 Nbnxm::InteractionLocality ilocality,
376 /*! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU */
377 void NbnxnDispatchKernel(nonbonded_verlet_t *nbv,
378 Nbnxm::InteractionLocality iLocality,
379 const interaction_const_t &ic,
383 gmx_enerdata_t *enerd,
386 #endif // GMX_NBNXN_NBNXN_H