2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 /*! \libinternal \file
38 * Functionality for per-atom data in the nbnxm module
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_nbnxm
46 #ifndef GMX_NBNXN_ATOMDATA_H
47 #define GMX_NBNXN_ATOMDATA_H
51 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
52 #include "gromacs/gpu_utils/hostallocator.h"
53 #include "gromacs/math/vectypes.h"
54 #include "gromacs/mdtypes/locality.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/bitmask.h"
57 #include "gromacs/utility/real.h"
65 struct nbnxn_atomdata_t;
66 struct nonbonded_verlet_t;
69 class GpuEventSynchronizer;
74 enum class KernelType;
77 //! Convenience type for vector with aligned memory
79 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
89 //! Stride for coordinate/force arrays with xyz coordinate storage
90 static constexpr int STRIDE_XYZ = 3;
91 //! Stride for coordinate/force arrays with xyzq coordinate storage
92 static constexpr int STRIDE_XYZQ = 4;
93 //! Size of packs of x, y or z with SIMD 4-grouped packed coordinates/forces
94 static constexpr int c_packX4 = 4;
95 //! Size of packs of x, y or z with SIMD 8-grouped packed coordinates/forces
96 static constexpr int c_packX8 = 8;
97 //! Stridefor a pack of 4 coordinates/forces
98 static constexpr int STRIDE_P4 = DIM * c_packX4;
99 //! Stridefor a pack of 8 coordinates/forces
100 static constexpr int STRIDE_P8 = DIM * c_packX8;
102 //! Returns the index in a coordinate array corresponding to atom a
103 template<int packSize>
104 static inline int atom_to_x_index(int a)
106 return DIM * (a & ~(packSize - 1)) + (a & (packSize - 1));
110 * \brief Struct that holds force and energy output buffers */
111 struct nbnxn_atomdata_output_t
113 /*! \brief Constructor
115 * \param[in] kernelType Type of non-bonded kernel
116 * \param[in] numEnergyGroups The number of energy groups
117 * \param[in] simdEnergyBufferStride Stride for entries in the energy buffers for SIMD kernels
118 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
120 nbnxn_atomdata_output_t(Nbnxm::KernelType kernelType,
122 int simdEnergyBufferStride,
123 gmx::PinningPolicy pinningPolicy);
125 //! f, size natoms*fstride
126 gmx::HostVector<real> f;
127 //! Shift force array, size SHIFTS*DIM
128 gmx::HostVector<real> fshift;
129 //! Temporary Van der Waals group energy storage
130 gmx::HostVector<real> Vvdw;
131 //! Temporary Coulomb group energy storage
132 gmx::HostVector<real> Vc;
133 //! Temporary SIMD Van der Waals group energy storage
134 AlignedVector<real> VSvdw;
135 //! Temporary SIMD Coulomb group energy storage
136 AlignedVector<real> VSc;
139 /*! \brief Block size in atoms for the non-bonded thread force-buffer reduction.
141 * Should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
142 * Should be small to reduce the reduction and zeroing cost,
143 * but too small will result in overhead.
144 * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
147 # define NBNXN_BUFFERFLAG_SIZE 8
149 # define NBNXN_BUFFERFLAG_SIZE 16
152 /*! \brief We store the reduction flags as gmx_bitmask_t.
153 * This limits the number of flags to BITMASK_SIZE.
155 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
158 //! LJ combination rules
159 enum class LJCombinationRule : int
163 //! Lorentz-Berthelot
171 //! String corresponding to LJ combination rule
172 const char* enumValueToString(LJCombinationRule enumValue);
175 * \brief Struct that stores atom related data for the nbnxn module
177 * Note: performance would improve slightly when all std::vector containers
178 * in this struct would not initialize during resize().
180 struct nbnxn_atomdata_t
181 { //NOLINT(clang-analyzer-optin.performance.Padding)
183 * \brief The actual atom data parameter values */
186 /*! \brief Constructor
188 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
190 Params(gmx::PinningPolicy pinningPolicy);
192 //! The number of different atom types
194 //! Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
195 gmx::HostVector<real> nbfp;
196 //! Combination rule, see enum defined above
197 LJCombinationRule ljCombinationRule;
198 //! LJ parameters per atom type, size numTypes*2
199 gmx::HostVector<real> nbfp_comb;
200 //! As nbfp, but with a stride for the present SIMD architecture
201 AlignedVector<real> nbfp_aligned;
202 //! Atom types per atom
203 gmx::HostVector<int> type;
204 //! LJ parameters per atom for fast SIMD loading
205 gmx::HostVector<real> lj_comb;
206 //! Charges per atom, not set with format nbatXYZQ
207 gmx::HostVector<real> q;
208 //! The number of energy groups
212 //! The energy groups, one int entry per cluster, only set when needed
213 gmx::HostVector<int> energrp;
217 * \brief Diagonal and topology exclusion helper data for all SIMD kernels. */
222 //! Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
223 AlignedVector<real> diagonal_4xn_j_minus_i;
224 //! Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
225 AlignedVector<real> diagonal_2xnn_j_minus_i;
226 //! Filters for topology exclusion masks for the SIMD kernels
227 AlignedVector<uint32_t> exclusion_filter;
228 //! Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
229 AlignedVector<uint64_t> exclusion_filter64;
230 //! Array of masks needed for exclusions
231 AlignedVector<real> interaction_array;
234 /*! \brief Constructor
236 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
238 nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
240 //! Returns a const reference to the parameters
241 const Params& params() const { return params_; }
243 //! Returns a non-const reference to the parameters
244 Params& paramsDeprecated() { return params_; }
246 //! Returns the current total number of atoms stored
247 int numAtoms() const { return numAtoms_; }
249 //! Return the coordinate buffer, and q with xFormat==nbatXYZQ
250 gmx::ArrayRef<const real> x() const { return x_; }
252 //! Return the coordinate buffer, and q with xFormat==nbatXYZQ
253 gmx::ArrayRef<real> x() { return x_; }
255 //! Resizes the coordinate buffer and sets the number of atoms
256 void resizeCoordinateBuffer(int numAtoms);
258 //! Resizes the force buffers for the current number of atoms
259 void resizeForceBuffers();
262 //! The LJ and charge parameters
264 //! The total number of atoms currently stored
268 //! Number of local atoms
270 //! The format of x (and q), enum
272 //! The format of f, enum
274 //! Do we need to update shift_vec every step?
275 gmx_bool bDynamicBox;
276 //! Shift vectors, copied from t_forcerec
277 gmx::HostVector<gmx::RVec> shift_vec;
278 //! stride for a coordinate in x (usually 3 or 4)
280 //! stride for a coordinate in f (usually 3 or 4)
284 //! x and possibly q, size natoms*xstride
285 gmx::HostVector<real> x_;
288 //! Masks for handling exclusions in the SIMD kernels
289 const SimdMasks simdMasks;
291 //! Output data structures, 1 per thread
292 std::vector<nbnxn_atomdata_output_t> out;
294 //! Reduction related data
296 //! Use the flags or operate on all atoms
297 gmx_bool bUseBufferFlags;
298 //! Flags for buffer zeroing+reduc.
299 std::vector<gmx_bitmask_t> buffer_flags;
300 //! Use tree for force reduction
301 gmx_bool bUseTreeReduce;
302 //! Synchronization step for tree reduce
303 tMPI_Atomic* syncStep;
307 /*! \brief Copy na rvec elements from x to xnb using nbatFormat, start dest a0,
308 * and fills up to na_round with coordinates that are far away.
310 void copy_rvec_to_nbat_real(const int* a, int na, int na_round, const rvec* x, int nbatFormat, real* xnb, int a0);
312 //! Describes the combination rule in use by this force field
315 enbnxninitcombruleDETECT,
316 enbnxninitcombruleGEOM,
317 enbnxninitcombruleLB,
318 enbnxninitcombruleNONE
321 /*! \brief Initialize the non-bonded atom data structure.
323 * The enum for nbatXFormat is in the file defining nbnxn_atomdata_t.
324 * Copy the ntypes*ntypes*2 sized nbfp non-bonded parameter list
325 * to the atom data structure.
326 * enbnxninitcombrule sets what combination rule data gets stored in nbat.
328 void nbnxn_atomdata_init(const gmx::MDLogger& mdlog,
329 nbnxn_atomdata_t* nbat,
330 Nbnxm::KernelType kernelType,
331 int enbnxninitcombrule,
333 gmx::ArrayRef<const real> nbfp,
337 //! Sets the atomdata after pair search
338 void nbnxn_atomdata_set(nbnxn_atomdata_t* nbat,
339 const Nbnxm::GridSet& gridSet,
340 gmx::ArrayRef<const int> atomTypes,
341 gmx::ArrayRef<const real> atomCharges,
342 gmx::ArrayRef<const int> atomInfo);
344 //! Copy the shift vectors to nbat
345 void nbnxn_atomdata_copy_shiftvec(gmx_bool dynamic_box, rvec* shift_vec, nbnxn_atomdata_t* nbat);
347 /*! \brief Transform coordinates to xbat layout
349 * Creates a copy of the coordinates buffer using short-range ordering.
351 * \param[in] gridSet The grids data.
352 * \param[in] locality If the transformation should be applied to local or non local coordinates.
353 * \param[in] fillLocal Tells if the local filler particle coordinates should be zeroed.
354 * \param[in] coordinates Coordinates in plain rvec format.
355 * \param[in,out] nbat Data in NBNXM format, used for mapping formats and to locate the output buffer.
357 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet& gridSet,
358 gmx::AtomLocality locality,
360 const rvec* coordinates,
361 nbnxn_atomdata_t* nbat);
363 /*! \brief Transform coordinates to xbat layout on GPU
365 * Creates a GPU copy of the coordinates buffer using short-range ordering.
366 * As input, uses coordinates in plain rvec format in GPU memory.
368 * \param[in] gridSet The grids data.
369 * \param[in] locality If the transformation should be applied to local or non local coordinates.
370 * \param[in] fillLocal Tells if the local filler particle coordinates should be zeroed.
371 * \param[in,out] gpu_nbv The NBNXM GPU data structure.
372 * \param[in] d_x Coordinates to be copied (in plain rvec format).
373 * \param[in] xReadyOnDevice Event synchronizer indicating that the coordinates are ready in the device memory.
375 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet& gridSet,
376 gmx::AtomLocality locality,
379 DeviceBuffer<gmx::RVec> d_x,
380 GpuEventSynchronizer* xReadyOnDevice);
382 /*! \brief Add the computed forces to \p f, an internal reduction might be performed as well
384 * \param[in] nbat Atom data in NBNXM format.
385 * \param[in] locality If the reduction should be performed on local or non-local atoms.
386 * \param[in] gridSet The grids data.
387 * \param[out] totalForce Buffer to accumulate resulting force
389 void reduceForces(nbnxn_atomdata_t* nbat, gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* totalForce);
391 //! Add the fshift force stored in nbat to fshift
392 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift);
394 //! Get the atom start index and number of atoms for a given locality
395 void nbnxn_get_atom_range(gmx::AtomLocality atomLocality,
396 const Nbnxm::GridSet& gridSet,