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/gpu_utils/devicebuffer_datatype.h"
104 #include "gromacs/math/vectypes.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 gmx_wallcycle;
120 struct interaction_const_t;
121 struct nonbonded_verlet_t;
132 /*! \brief Switch for whether to use GPU for buffer ops*/
133 enum class BufferOpsUseGpu
139 class GpuEventSynchronizer;
143 class ForceWithShiftForces;
145 class UpdateGroupsCog;
150 enum class KernelType;
156 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
157 enum class KernelType : int
168 /*! \brief Ewald exclusion types */
169 enum class EwaldExclusionType : int
177 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
180 //! The non-bonded type, also affects the pairlist construction kernel
181 KernelType kernelType = KernelType::NotSet;
182 //! Ewald exclusion computation handling type, currently only used for CPU
183 EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
186 /*! \brief Return a string identifying the kernel type.
188 * \param [in] kernelType nonbonded kernel type, takes values from the nbnxn_kernel_type enum
189 * \returns a string identifying the kernel corresponding to the type passed as argument
191 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
195 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
197 enbvClearFNo, enbvClearFYes
201 * \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
202 struct nonbonded_verlet_t
205 //! Constructs an object from its components
206 nonbonded_verlet_t(std::unique_ptr<PairlistSets> pairlistSets,
207 std::unique_ptr<PairSearch> pairSearch,
208 std::unique_ptr<nbnxn_atomdata_t> nbat,
209 const Nbnxm::KernelSetup &kernelSetup,
210 gmx_nbnxn_gpu_t *gpu_nbv,
211 gmx_wallcycle *wcycle);
213 ~nonbonded_verlet_t();
215 //! Returns whether a GPU is use for the non-bonded calculations
218 return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
221 //! Returns whether a GPU is emulated for the non-bonded calculations
222 bool emulateGpu() const
224 return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
227 //! Return whether the pairlist is of simple, CPU type
228 bool pairlistIsSimple() const
230 return !useGpu() && !emulateGpu();
233 //! Initialize the pair list sets, TODO this should be private
234 void initPairlistSets(bool haveMultipleDomains);
236 //! Returns the order of the local atoms on the grid
237 gmx::ArrayRef<const int> getLocalAtomOrder() const;
239 //! Sets the order of the local atoms to the order grid atom ordering
240 void setLocalAtomOrder();
242 //! Returns the index position of the atoms on the search grid
243 gmx::ArrayRef<const int> getGridIndices() const;
245 //! Constructs the pairlist for the given locality
246 void constructPairlist(Nbnxm::InteractionLocality iLocality,
247 const t_blocka *excl,
251 //! Updates all the atom properties in Nbnxm
252 void setAtomProperties(const t_mdatoms &mdatoms,
253 gmx::ArrayRef<const int> atomInfo);
255 /*!\brief Convert the coordinates to NBNXM format for the given locality.
257 * The API function for the transformation of the coordinates from one layout to another.
259 * \param[in] locality Whether coordinates for local or non-local atoms should be transformed.
260 * \param[in] fillLocal If the coordinates for filler particles should be zeroed.
261 * \param[in] coordinates Coordinates in plain rvec format to be transformed.
263 void convertCoordinates(Nbnxm::AtomLocality locality,
265 gmx::ArrayRef<const gmx::RVec> coordinates);
267 /*!\brief Convert the coordinates to NBNXM format on the GPU for the given locality
269 * The API function for the transformation of the coordinates from one layout to another in the GPU memory.
271 * \param[in] locality Whether coordinates for local or non-local atoms should be transformed.
272 * \param[in] fillLocal If the coordinates for filler particles should be zeroed.
273 * \param[in] d_x GPU coordinates buffer in plain rvec format to be transformed.
275 void convertCoordinatesGpu(Nbnxm::AtomLocality locality,
277 DeviceBuffer<float> d_x);
279 //! Init for GPU version of setup coordinates in Nbnxm
280 void atomdata_init_copy_x_to_nbat_x_gpu();
282 //! Sync the nonlocal GPU stream with dependent tasks in the local queue.
283 void insertNonlocalGpuDependency(Nbnxm::InteractionLocality interactionLocality);
285 //! Returns a reference to the pairlist sets
286 const PairlistSets &pairlistSets() const
288 return *pairlistSets_;
291 //! Returns whether step is a dynamic list pruning step, for CPU lists
292 bool isDynamicPruningStepCpu(int64_t step) const;
294 //! Returns whether step is a dynamic list pruning step, for GPU lists
295 bool isDynamicPruningStepGpu(int64_t step) const;
297 //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
298 void dispatchPruneKernelCpu(Nbnxm::InteractionLocality iLocality,
299 const rvec *shift_vec);
301 //! Dispatches the dynamic pruning kernel for GPU lists
302 void dispatchPruneKernelGpu(int64_t step);
304 //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
305 void dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality,
306 const interaction_const_t &ic,
307 const gmx::StepWorkload &stepWork,
309 const t_forcerec &fr,
310 gmx_enerdata_t *enerd,
313 //! Executes the non-bonded free-energy kernel, always runs on the CPU
314 void dispatchFreeEnergyKernel(Nbnxm::InteractionLocality iLocality,
315 const t_forcerec *fr,
317 gmx::ForceWithShiftForces *forceWithShiftForces,
318 const t_mdatoms &mdatoms,
321 gmx_enerdata_t *enerd,
322 const gmx::StepWorkload &stepWork,
325 /*! \brief Add the forces stored in nbat to f, zeros the forces in nbat
326 * \param [in] locality Local or non-local
327 * \param [inout] force Force to be added to
329 void atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality locality,
330 gmx::ArrayRef<gmx::RVec> force);
332 /*! \brief Add the forces stored in nbat to total force using GPU buffer opse
334 * \param [in] locality Local or non-local
335 * \param [in,out] totalForcesDevice Force to be added to
336 * \param [in] forcesPmeDevice Device buffer with PME forces
337 * \param [in] pmeForcesReady Event triggered when PME force calculation has completed
338 * \param [in] useGpuFPmeReduction Whether PME forces should be added
339 * \param [in] accumulateForce If the total force buffer already contains data
341 void atomdata_add_nbat_f_to_f_gpu(Nbnxm::AtomLocality locality,
342 DeviceBuffer<float> totalForcesDevice,
343 void *forcesPmeDevice,
344 GpuEventSynchronizer *pmeForcesReady,
345 bool useGpuFPmeReduction,
346 bool accumulateForce);
348 /*! \brief Outer body of function to perform initialization for F buffer operations on GPU. */
349 void atomdata_init_add_nbat_f_to_f_gpu();
351 /*! \brief Wait for GPU force reduction task and D2H transfer of its results to complete
353 * FIXME: need more details: when should be called / after which operation, etc.
355 void wait_for_gpu_force_reduction(Nbnxm::AtomLocality locality);
357 /*! \brief return pointer to GPU event recorded when coordinates have been copied to device */
358 void* get_x_on_device_event();
360 /*! \brief Wait for non-local copy of coordinate buffer from device to host */
361 void wait_nonlocal_x_copy_D2H_done();
363 /*! \brief return GPU pointer to f in rvec format */
364 void* get_gpu_frvec();
366 /*! \brief Ensure local stream waits for non-local stream */
367 void stream_local_wait_for_nonlocal();
369 //! Return the kernel setup
370 const Nbnxm::KernelSetup &kernelSetup() const
375 //! Returns the outer radius for the pair list
376 real pairlistInnerRadius() const;
378 //! Returns the outer radius for the pair list
379 real pairlistOuterRadius() const;
381 //! Changes the pair-list outer and inner radius
382 void changePairlistRadii(real rlistOuter,
385 //! Set up internal flags that indicate what type of short-range work there is.
386 void setupGpuShortRangeWork(const gmx::GpuBonded *gpuBonded,
387 const Nbnxm::InteractionLocality iLocality)
389 if (useGpu() && !emulateGpu())
391 Nbnxm::setupGpuShortRangeWork(gpu_nbv, gpuBonded, iLocality);
395 //! Returns true if there is GPU short-range work for the given atom locality.
396 bool haveGpuShortRangeWork(const Nbnxm::AtomLocality aLocality)
398 return ((useGpu() && !emulateGpu()) &&
399 Nbnxm::haveGpuShortRangeWork(gpu_nbv, aLocality));
402 // TODO: Make all data members private
404 //! All data related to the pair lists
405 std::unique_ptr<PairlistSets> pairlistSets_;
406 //! Working data for constructing the pairlists
407 std::unique_ptr<PairSearch> pairSearch_;
409 std::unique_ptr<nbnxn_atomdata_t> nbat;
411 //! The non-bonded setup, also affects the pairlist construction kernel
412 Nbnxm::KernelSetup kernelSetup_;
413 //! \brief Pointer to wallcycle structure.
414 gmx_wallcycle *wcycle_;
416 //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
417 gmx_nbnxn_gpu_t *gpu_nbv;
423 /*! \brief Creates an Nbnxm object */
424 std::unique_ptr<nonbonded_verlet_t>
425 init_nb_verlet(const gmx::MDLogger &mdlog,
426 gmx_bool bFEP_NonBonded,
427 const t_inputrec *ir,
428 const t_forcerec *fr,
430 const gmx_hw_info_t &hardwareInfo,
431 const gmx_device_info_t *deviceInfo,
432 const gmx_mtop_t *mtop,
434 gmx_wallcycle *wcycle);
438 /*! \brief Put the atoms on the pair search grid.
440 * Only atoms atomStart to atomEnd in x are put on the grid.
441 * The atom_density is used to determine the grid size.
442 * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
443 * With domain decomposition part of the n particles might have migrated,
444 * but have not been removed yet. This count is given by nmoved.
445 * When move[i] < 0 particle i has migrated and will not be put on the grid.
446 * Without domain decomposition move will be NULL.
448 void nbnxn_put_on_grid(nonbonded_verlet_t *nb_verlet,
451 const rvec lowerCorner,
452 const rvec upperCorner,
453 const gmx::UpdateGroupsCog *updateGroupsCog,
457 gmx::ArrayRef<const int> atomInfo,
458 gmx::ArrayRef<const gmx::RVec> x,
462 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
464 * with domain decomposition. Should be called after calling
465 * nbnxn_search_put_on_grid for the local atoms / home zone.
467 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t *nb_verlet,
468 const struct gmx_domdec_zones_t *zones,
469 gmx::ArrayRef<const int> atomInfo,
470 gmx::ArrayRef<const gmx::RVec> x);
472 #endif // GMX_NBNXN_NBNXM_H