2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,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 #ifndef GMX_NBNXN_ATOMDATA_H
37 #define GMX_NBNXN_ATOMDATA_H
41 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
42 #include "gromacs/gpu_utils/hostallocator.h"
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdtypes/locality.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/bitmask.h"
47 #include "gromacs/utility/real.h"
49 #include "gpu_types.h"
56 struct nbnxn_atomdata_t;
57 struct nonbonded_verlet_t;
61 class GpuEventSynchronizer;
66 enum class KernelType;
69 /* Convenience type for vector with aligned memory */
71 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
74 nbatXYZ, nbatXYZQ, nbatX4, nbatX8
77 //! Stride for coordinate/force arrays with xyz coordinate storage
78 static constexpr int STRIDE_XYZ = 3;
79 //! Stride for coordinate/force arrays with xyzq coordinate storage
80 static constexpr int STRIDE_XYZQ = 4;
81 //! Size of packs of x, y or z with SIMD 4-grouped packed coordinates/forces
82 static constexpr int c_packX4 = 4;
83 //! Size of packs of x, y or z with SIMD 8-grouped packed coordinates/forces
84 static constexpr int c_packX8 = 8;
85 //! Stridefor a pack of 4 coordinates/forces
86 static constexpr int STRIDE_P4 = DIM*c_packX4;
87 //! Stridefor a pack of 8 coordinates/forces
88 static constexpr int STRIDE_P8 = DIM*c_packX8;
90 //! Returns the index in a coordinate array corresponding to atom a
91 template<int packSize> static inline int atom_to_x_index(int a)
93 return DIM*(a & ~(packSize - 1)) + (a & (packSize - 1));
96 // Struct that holds force and energy output buffers
97 struct nbnxn_atomdata_output_t
101 * \param[in] kernelType Type of non-bonded kernel
102 * \param[in] numEnergyGroups The number of energy groups
103 * \param[in] simdEnergyBufferStride Stride for entries in the energy buffers for SIMD kernels
104 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
106 nbnxn_atomdata_output_t(Nbnxm::KernelType kernelType,
108 int simdEnergyBUfferStride,
109 gmx::PinningPolicy pinningPolicy);
111 gmx::HostVector<real> f; // f, size natoms*fstride
112 gmx::HostVector<real> fshift; // Shift force array, size SHIFTS*DIM
113 gmx::HostVector<real> Vvdw; // Temporary Van der Waals group energy storage
114 gmx::HostVector<real> Vc; // Temporary Coulomb group energy storage
115 AlignedVector<real> VSvdw; // Temporary SIMD Van der Waals group energy storage
116 AlignedVector<real> VSc; // Temporary SIMD Coulomb group energy storage
119 /* Block size in atoms for the non-bonded thread force-buffer reduction,
120 * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
121 * Should be small to reduce the reduction and zeroing cost,
122 * but too small will result in overhead.
123 * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
126 #define NBNXN_BUFFERFLAG_SIZE 8
128 #define NBNXN_BUFFERFLAG_SIZE 16
131 /* We store the reduction flags as gmx_bitmask_t.
132 * This limits the number of flags to BITMASK_SIZE.
134 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
136 /* Flags for telling if threads write to force output buffers */
138 int nflag; /* The number of flag blocks */
139 gmx_bitmask_t *flag; /* Bit i is set when thread i writes to a cell-block */
140 int flag_nalloc; /* Allocation size of cxy_flag */
141 } nbnxn_buffer_flags_t;
143 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
145 ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
148 /* Struct that stores atom related data for the nbnxn module
150 * Note: performance would improve slightly when all std::vector containers
151 * in this struct would not initialize during resize().
153 struct nbnxn_atomdata_t
154 { //NOLINT(clang-analyzer-optin.performance.Padding)
159 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
161 Params(gmx::PinningPolicy pinningPolicy);
163 // The number of different atom types
165 // Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
166 gmx::HostVector<real> nbfp;
167 // Combination rule, see enum defined above
169 // LJ parameters per atom type, size numTypes*2
170 gmx::HostVector<real> nbfp_comb;
171 // As nbfp, but with a stride for the present SIMD architecture
172 AlignedVector<real> nbfp_aligned;
173 // Atom types per atom
174 gmx::HostVector<int> type;
175 // LJ parameters per atom for fast SIMD loading
176 gmx::HostVector<real> lj_comb;
177 // Charges per atom, not set with format nbatXYZQ
178 gmx::HostVector<real> q;
179 // The number of energy groups
183 // The energy groups, one int entry per cluster, only set when needed
184 gmx::HostVector<int> energrp;
187 // Diagonal and topology exclusion helper data for all SIMD kernels
192 // Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
193 AlignedVector<real> diagonal_4xn_j_minus_i;
194 // Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
195 AlignedVector<real> diagonal_2xnn_j_minus_i;
196 // Filters for topology exclusion masks for the SIMD kernels
197 AlignedVector<uint32_t> exclusion_filter;
198 // Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
199 AlignedVector<uint64_t> exclusion_filter64;
200 // Array of masks needed for exclusions
201 AlignedVector<real> interaction_array;
206 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
208 nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
210 /* Returns a const reference to the parameters */
211 const Params ¶ms() const
216 /* Returns a non-const reference to the parameters */
217 Params ¶msDeprecated()
222 /* Returns the current total number of atoms stored */
228 /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
229 gmx::ArrayRef<const real> x() const
234 /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
235 gmx::ArrayRef<real> x()
240 /* Resizes the coordinate buffer and sets the number of atoms */
241 void resizeCoordinateBuffer(int numAtoms);
243 /* Resizes the force buffers for the current number of atoms */
244 void resizeForceBuffers();
247 // The LJ and charge parameters
249 // The total number of atoms currently stored
252 int natoms_local; /* Number of local atoms */
253 int XFormat; /* The format of x (and q), enum */
254 int FFormat; /* The format of f, enum */
255 gmx_bool bDynamicBox; /* Do we need to update shift_vec every step? */
256 gmx::HostVector<gmx::RVec> shift_vec; /* Shift vectors, copied from t_forcerec */
257 int xstride; /* stride for a coordinate in x (usually 3 or 4) */
258 int fstride; /* stride for a coordinate in f (usually 3 or 4) */
260 gmx::HostVector<real> x_; /* x and possibly q, size natoms*xstride */
263 // Masks for handling exclusions in the SIMD kernels
264 const SimdMasks simdMasks;
267 std::vector<nbnxn_atomdata_output_t> out; /* Output data structures, 1 per thread */
269 /* Reduction related data */
270 gmx_bool bUseBufferFlags; /* Use the flags or operate on all atoms */
271 nbnxn_buffer_flags_t buffer_flags; /* Flags for buffer zeroing+reduc. */
272 gmx_bool bUseTreeReduce; /* Use tree for force reduction */
273 tMPI_Atomic *syncStep; /* Synchronization step for tree reduce */
276 /* Copy na rvec elements from x to xnb using nbatFormat, start dest a0,
277 * and fills up to na_round with coordinates that are far away.
279 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
280 const rvec *x, int nbatFormat,
284 enbnxninitcombruleDETECT, enbnxninitcombruleGEOM, enbnxninitcombruleLB, enbnxninitcombruleNONE
287 /* Initialize the non-bonded atom data structure.
288 * The enum for nbatXFormat is in the file defining nbnxn_atomdata_t.
289 * Copy the ntypes*ntypes*2 sized nbfp non-bonded parameter list
290 * to the atom data structure.
291 * enbnxninitcombrule sets what combination rule data gets stored in nbat.
293 void nbnxn_atomdata_init(const gmx::MDLogger &mdlog,
294 nbnxn_atomdata_t *nbat,
295 Nbnxm::KernelType kernelType,
296 int enbnxninitcombrule,
297 int ntype, const real *nbfp,
301 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
302 const Nbnxm::GridSet &gridSet,
303 const t_mdatoms *mdatoms,
306 /* Copy the shift vectors to nbat */
307 void nbnxn_atomdata_copy_shiftvec(gmx_bool dynamic_box,
309 nbnxn_atomdata_t *nbat);
311 /*! \brief Transform coordinates to xbat layout
313 * Creates a copy of the coordinates buffer using short-range ordering.
315 * \param[in] gridSet The grids data.
316 * \param[in] locality If the transformation should be applied to local or non local coordinates.
317 * \param[in] fillLocal Tells if the local filler particle coordinates should be zeroed.
318 * \param[in] coordinates Coordinates in plain rvec format.
319 * \param[in,out] nbat Data in NBNXM format, used for mapping formats and to locate the output buffer.
321 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet &gridSet,
322 gmx::AtomLocality locality,
324 const rvec *coordinates,
325 nbnxn_atomdata_t *nbat);
327 /*! \brief Transform coordinates to xbat layout on GPU
329 * Creates a GPU copy of the coordinates buffer using short-range ordering.
330 * As input, uses coordinates in plain rvec format in GPU memory.
332 * \param[in] gridSet The grids data.
333 * \param[in] locality If the transformation should be applied to local or non local coordinates.
334 * \param[in] fillLocal Tells if the local filler particle coordinates should be zeroed.
335 * \param[in,out] gpu_nbv The NBNXM GPU data structure.
336 * \param[in] d_x Coordinates to be copied (in plain rvec format).
337 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
339 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet &gridSet,
340 gmx::AtomLocality locality,
342 gmx_nbnxn_gpu_t *gpu_nbv,
343 DeviceBuffer<float> d_x,
344 GpuEventSynchronizer *xReadyOnDevice);
346 /*! \brief Add the computed forces to \p f, an internal reduction might be performed as well
348 * \param[in] nbat Atom data in NBNXM format.
349 * \param[in] locality If the reduction should be performed on local or non-local atoms.
350 * \param[in] gridSet The grids data.
351 * \param[out] totalForce Buffer to accumulate resulting force
353 void reduceForces(nbnxn_atomdata_t *nbat,
354 gmx::AtomLocality locality,
355 const Nbnxm::GridSet &gridSet,
358 /*! \brief Reduce forces on the GPU
360 * \param[in] locality If the reduction should be performed on local or non-local atoms.
361 * \param[out] totalForcesDevice Device buffer to accumulate resulting force.
362 * \param[in] gridSet The grids data.
363 * \param[in] pmeForcesDevice Device buffer with PME forces.
364 * \param[in] dependencyList List of synchronizers that represent the dependencies the reduction task needs to sync on.
365 * \param[in] gpu_nbv The NBNXM GPU data structure.
366 * \param[in] useGpuFPmeReduction Whether PME forces should be added.
367 * \param[in] accumulateForce Whether there are usefull data already in the total force buffer.
369 void reduceForcesGpu(gmx::AtomLocality locality,
370 DeviceBuffer<float> totalForcesDevice,
371 const Nbnxm::GridSet &gridSet,
372 void *pmeForcesDevice,
373 gmx::ArrayRef<GpuEventSynchronizer* const> dependencyList,
374 gmx_nbnxn_gpu_t *gpu_nbv,
375 bool useGpuFPmeReduction,
376 bool accumulateForce);
378 /* Add the fshift force stored in nbat to fshift */
379 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t &nbat,
380 gmx::ArrayRef<gmx::RVec> fshift);
382 /* Get the atom start index and number of atoms for a given locality */
383 void nbnxn_get_atom_range(gmx::AtomLocality atomLocality,
384 const Nbnxm::GridSet &gridSet,