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/utility/arrayref.h"
105 #include "gromacs/utility/enumerationhelpers.h"
106 #include "gromacs/utility/real.h"
108 #include "locality.h"
110 // TODO: Remove this include and the two nbnxm includes above
111 #include "nbnxm_gpu.h"
113 struct gmx_device_info_t;
114 struct gmx_domdec_zones_t;
115 struct gmx_enerdata_t;
116 struct gmx_hw_info_t;
118 struct gmx_wallcycle;
119 struct interaction_const_t;
120 struct nbnxn_pairlist_set_t;
122 struct nonbonded_verlet_t;
123 enum class PairlistType;
135 class UpdateGroupsCog;
140 enum class KernelType;
144 * \brief The setup for generating and pruning the nbnxn pair list.
146 * Without dynamic pruning rlistOuter=rlistInner.
148 struct NbnxnListParameters
150 /*! \brief Constructor producing a struct with dynamic pruning disabled
152 NbnxnListParameters(Nbnxm::KernelType kernelType,
154 bool haveMultipleDomains);
156 PairlistType pairlistType; //!< The type of cluster-pair list
157 real rlistOuter; //!< Cut-off of the larger, outer pair-list
158 real rlistInner; //!< Cut-off of the smaller, inner pair-list
159 bool haveMultipleDomains; //!< True when using DD with multiple domains
160 bool useDynamicPruning; //!< Are we using dynamic pair-list pruning
161 int nstlistPrune; //!< Pair-list dynamic pruning interval
162 int numRollingPruningParts; //!< The number parts to divide the pair-list into for rolling pruning, a value of 1 gives no rolling pruning
163 int lifetime; //!< Lifetime in steps of the pair-list
166 /*! \brief Resources that can be used to execute non-bonded kernels on */
167 enum class NonbondedResource : int
177 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
178 enum class KernelType : int
189 /*! \brief Ewald exclusion types */
190 enum class EwaldExclusionType : int
198 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
201 //! The non-bonded type, also affects the pairlist construction kernel
202 KernelType kernelType = KernelType::NotSet;
203 //! Ewald exclusion computation handling type, currently only used for CPU
204 EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
207 /*! \brief Return a string identifying the kernel type.
209 * \param [in] kernelType nonbonded kernel type, takes values from the nbnxn_kernel_type enum
210 * \returns a string identifying the kernel corresponding to the type passed as argument
212 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
216 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
218 enbvClearFNo, enbvClearFYes
221 /*! \brief Generates a pair-list for the given locality.
223 * With perturbed particles, also a group scheme style nbl_fep list is made.
225 void nbnxn_make_pairlist(nonbonded_verlet_t *nbv,
226 Nbnxm::InteractionLocality iLocality,
227 nbnxn_pairlist_set_t *pairlistSet,
228 const t_blocka *excl,
232 /*! \brief Prune all pair-lists with given locality (currently CPU only)
234 * For all pair-lists with given locality, takes the outer list and prunes out
235 * pairs beyond the pairlist inner radius and writes the result to a list that is
236 * to be consumed by the non-bonded kernel.
238 void NbnxnDispatchPruneKernel(nbnxn_pairlist_set_t *pairlistSet,
239 Nbnxm::KernelType kernelType,
240 const nbnxn_atomdata_t *nbat,
241 const rvec *shift_vec);
244 * \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
245 struct nonbonded_verlet_t
251 PairlistSets(const NbnxnListParameters &listParams,
252 bool haveMultipleDomains,
253 int minimumIlistCountForGpuBalancing);
255 //! Construct the pairlist set for the given locality
256 void construct(Nbnxm::InteractionLocality iLocality,
258 nbnxn_atomdata_t *nbat,
259 const t_blocka *excl,
260 Nbnxm::KernelType kernelbType,
264 //! Dispatches the dynamic pruning kernel for the given locality
265 void dispatchPruneKernel(Nbnxm::InteractionLocality iLocality,
266 const nbnxn_atomdata_t *nbat,
267 const rvec *shift_vec,
268 Nbnxm::KernelType kernelbType);
270 //! Returns the pair list parameters
271 const NbnxnListParameters ¶ms() const
276 //! Returns the number of steps performed with the current pair list
277 int numStepsWithPairlist(int64_t step) const
279 return step - outerListCreationStep_;
282 //! Returns whether step is a dynamic list pruning step, for CPU lists
283 bool isDynamicPruningStepCpu(int64_t step) const
285 return (params_.useDynamicPruning &&
286 numStepsWithPairlist(step) % params_.nstlistPrune == 0);
289 //! Returns whether step is a dynamic list pruning step, for GPU lists
290 bool isDynamicPruningStepGpu(int64_t step) const
292 const int age = numStepsWithPairlist(step);
294 return (params_.useDynamicPruning &&
296 age < params_.lifetime &&
297 (params_.haveMultipleDomains || age % 2 == 0));
300 //! Changes the pair-list outer and inner radius
301 void changeRadii(real rlistOuter,
304 params_.rlistOuter = rlistOuter;
305 params_.rlistInner = rlistInner;
308 //! Returns the pair-list set for the given locality
309 const nbnxn_pairlist_set_t &pairlistSet(Nbnxm::InteractionLocality iLocality) const
311 if (iLocality == Nbnxm::InteractionLocality::Local)
317 GMX_ASSERT(nonlocalSet_, "Need a non-local set when requesting access");
318 return *nonlocalSet_;
323 //! Returns the pair-list set for the given locality
324 nbnxn_pairlist_set_t &pairlistSet(Nbnxm::InteractionLocality iLocality)
326 if (iLocality == Nbnxm::InteractionLocality::Local)
332 GMX_ASSERT(nonlocalSet_, "Need a non-local set when requesting access");
333 return *nonlocalSet_;
337 //! Parameters for the search and list pruning setup
338 NbnxnListParameters params_;
339 //! Pair list balancing parameter for use with GPU
340 int minimumIlistCountForGpuBalancing_;
341 //! Local pairlist set
342 std::unique_ptr<nbnxn_pairlist_set_t> localSet_;
343 //! Non-local pairlist set
344 std::unique_ptr<nbnxn_pairlist_set_t> nonlocalSet_;
345 //! MD step at with the outer lists in pairlistSets_ were created
346 int64_t outerListCreationStep_;
349 //! Constructs an object from its components
350 nonbonded_verlet_t(std::unique_ptr<PairlistSets> pairlistSets_,
351 std::unique_ptr<nbnxn_search> nbs,
352 std::unique_ptr<nbnxn_atomdata_t> nbat,
353 const Nbnxm::KernelSetup &kernelSetup,
354 gmx_nbnxn_gpu_t *gpu_nbv);
356 ~nonbonded_verlet_t();
358 //! Returns whether a GPU is use for the non-bonded calculations
361 return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
364 //! Returns whether a GPU is emulated for the non-bonded calculations
365 bool emulateGpu() const
367 return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
370 //! Return whether the pairlist is of simple, CPU type
371 bool pairlistIsSimple() const
373 return !useGpu() && !emulateGpu();
376 //! Initialize the pair list sets, TODO this should be private
377 void initPairlistSets(bool haveMultipleDomains);
379 //! Sets the order of the local atoms to the order grid atom ordering
380 void setLocalAtomOrder();
382 //! Constructs the pairlist for the given locality
383 void constructPairlist(Nbnxm::InteractionLocality iLocality,
384 const t_blocka *excl,
388 //! Returns a reference to the pairlist sets
389 const PairlistSets &pairlistSets() const
391 return *pairlistSets_;
394 //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
395 void dispatchPruneKernelCpu(Nbnxm::InteractionLocality iLocality,
396 const rvec *shift_vec);
398 //! Dispatches the dynamic pruning kernel for GPU lists
399 void dispatchPruneKernelGpu(int64_t step)
401 const bool stepIsEven = (pairlistSets().numStepsWithPairlist(step) % 2 == 0);
403 Nbnxm::gpu_launch_kernel_pruneonly(gpu_nbv,
404 stepIsEven ? Nbnxm::InteractionLocality::Local : Nbnxm::InteractionLocality::NonLocal,
405 pairlistSets().params().numRollingPruningParts);
408 //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
409 void dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality,
410 const interaction_const_t &ic,
414 gmx_enerdata_t *enerd,
417 //! Executes the non-bonded free-energy kernel, always runs on the CPU
418 void dispatchFreeEnergyKernel(Nbnxm::InteractionLocality iLocality,
422 const t_mdatoms &mdatoms,
425 gmx_enerdata_t *enerd,
429 //! Add the forces stored in nbat to f, zeros the forces in nbat */
430 void atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality locality,
432 gmx_wallcycle *wcycle);
434 //! Return the kernel setup
435 const Nbnxm::KernelSetup &kernelSetup() const
440 /*! \brief Returns the number of x and y cells in the local grid */
441 void getLocalNumCells(int *numCellsX,
442 int *numCellsY) const;
444 // TODO: Make all data members private
446 //! All data related to the pair lists
447 std::unique_ptr<PairlistSets> pairlistSets_;
448 //! Working data for constructing the pairlists
449 std::unique_ptr<nbnxn_search> nbs;
451 std::unique_ptr<nbnxn_atomdata_t> nbat;
453 //! The non-bonded setup, also affects the pairlist construction kernel
454 Nbnxm::KernelSetup kernelSetup_;
456 //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
457 gmx_nbnxn_gpu_t *gpu_nbv;
463 /*! \brief Creates an Nbnxm object */
464 std::unique_ptr<nonbonded_verlet_t>
465 init_nb_verlet(const gmx::MDLogger &mdlog,
466 gmx_bool bFEP_NonBonded,
467 const t_inputrec *ir,
468 const t_forcerec *fr,
470 const gmx_hw_info_t &hardwareInfo,
471 const gmx_device_info_t *deviceInfo,
472 const gmx_mtop_t *mtop,
477 /*! \brief Put the atoms on the pair search grid.
479 * Only atoms atomStart to atomEnd in x are put on the grid.
480 * The atom_density is used to determine the grid size.
481 * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
482 * With domain decomposition part of the n particles might have migrated,
483 * but have not been removed yet. This count is given by nmoved.
484 * When move[i] < 0 particle i has migrated and will not be put on the grid.
485 * Without domain decomposition move will be NULL.
487 void nbnxn_put_on_grid(nonbonded_verlet_t *nb_verlet,
490 const rvec lowerCorner,
491 const rvec upperCorner,
492 const gmx::UpdateGroupsCog *updateGroupsCog,
497 gmx::ArrayRef<const gmx::RVec> x,
501 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
503 * with domain decomposition. Should be called after calling
504 * nbnxn_search_put_on_grid for the local atoms / home zone.
506 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t *nb_verlet,
507 const struct gmx_domdec_zones_t *zones,
509 gmx::ArrayRef<const gmx::RVec> x);
511 /*! \brief Returns the order indices of the atoms on the pairlist search grid */
512 gmx::ArrayRef<const int> nbnxn_get_atomorder(const nbnxn_search* nbs);
514 /*! \brief Renumbers the atom indices on the grid to consecutive order */
515 void nbnxn_set_atomorder(nbnxn_search *nbs);
517 /*! \brief Returns the index position of the atoms on the pairlist search grid */
518 gmx::ArrayRef<const int> nbnxn_get_gridindices(const nbnxn_search* nbs);
520 #endif // GMX_NBNXN_NBNXN_H