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_NBNXM_PAIRLIST_H
37 #define GMX_NBNXM_PAIRLIST_H
43 #include "gromacs/gpu_utils/hostallocator.h"
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/mdtypes/nblist.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/bitmask.h"
48 #include "gromacs/utility/defaultinitializationallocator.h"
49 #include "gromacs/utility/enumerationhelpers.h"
50 #include "gromacs/utility/real.h"
52 // This file with constants is separate from this file to be able
53 // to include it during OpenCL jitting without including config.h
54 #include "gromacs/nbnxm/constants.h"
56 struct NbnxnPairlistCpuWork;
57 struct NbnxnPairlistGpuWork;
62 enum class KernelType;
65 /* Convenience type for vector with aligned memory */
67 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
69 /* Convenience type for vector that avoids initialization at resize() */
71 using FastVector = std::vector < T, gmx::DefaultInitializationAllocator < T>>;
73 enum class PairlistType : int
82 static constexpr gmx::EnumerationArray<PairlistType, int> IClusterSizePerListType = { 4, 4, 4, 8 };
83 static constexpr gmx::EnumerationArray<PairlistType, int> JClusterSizePerListType = { 2, 4, 8, 8 };
87 /*! \brief The setup for generating and pruning the nbnxn pair list.
89 * Without dynamic pruning rlistOuter=rlistInner.
91 struct NbnxnListParameters
93 /*! \brief Constructor producing a struct with dynamic pruning disabled
95 NbnxnListParameters(Nbnxm::KernelType kernelType,
98 PairlistType pairlistType; //!< The type of cluster-pair list
99 bool useDynamicPruning; //!< Are we using dynamic pair-list pruning
100 int nstlistPrune; //!< Pair-list dynamic pruning interval
101 real rlistOuter; //!< Cut-off of the larger, outer pair-list
102 real rlistInner; //!< Cut-off of the smaller, inner pair-list
103 int numRollingParts; //!< The number parts to divide the pair-list into for rolling pruning, a value of 1 gives no rolling pruning
108 /* With CPU kernels the i-cluster size is always 4 atoms. */
109 static constexpr int c_nbnxnCpuIClusterSize = 4;
111 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
112 #if GMX_GPU == GMX_GPU_OPENCL
113 static constexpr int c_nbnxnGpuClusterSize = GMX_OPENCL_NB_CLUSTER_SIZE;
115 static constexpr int c_nbnxnGpuClusterSize = 8;
118 /* The number of clusters in a pair-search cell, used for GPU */
119 static constexpr int c_gpuNumClusterPerCellZ = 2;
120 static constexpr int c_gpuNumClusterPerCellY = 2;
121 static constexpr int c_gpuNumClusterPerCellX = 2;
122 static constexpr int c_gpuNumClusterPerCell = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
125 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
126 * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
128 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
130 /* The fixed size of the exclusion mask array for a half cluster pair */
131 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
133 /* A buffer data structure of 64 bytes
134 * to be placed at the beginning and end of structs
135 * to avoid cache invalidation of the real contents
136 * of the struct by writes to neighboring memory.
140 } gmx_cache_protect_t;
142 /* Abstract type for pair searching data */
143 typedef struct nbnxn_search * nbnxn_search_t;
145 /* Function that should return a pointer *ptr to memory
147 * Error handling should be done within this function.
149 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
151 /* Function that should free the memory pointed to by *ptr.
152 * NULL should not be passed to this function.
154 typedef void nbnxn_free_t (void *ptr);
156 /* This is the actual cluster-pair list j-entry.
157 * cj is the j-cluster.
158 * The interaction bits in excl are indexed i-major, j-minor.
159 * The cj entries are sorted such that ones with exclusions come first.
160 * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
161 * is found, all subsequent j-entries in the i-entry also have full masks.
165 int cj; /* The j-cluster */
166 unsigned int excl; /* The exclusion (interaction) bits */
169 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
170 * The upper bits contain information for non-bonded kernel optimization.
171 * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
172 * But three flags can be used to skip interactions, currently only for subc=0
173 * !(shift & NBNXN_CI_DO_LJ(subc)) => we can skip LJ for all pairs
174 * shift & NBNXN_CI_HALF_LJ(subc) => we can skip LJ for the second half of i
175 * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
177 #define NBNXN_CI_SHIFT 127
178 #define NBNXN_CI_DO_LJ(subc) (1<<(7+3*(subc)))
179 #define NBNXN_CI_HALF_LJ(subc) (1<<(8+3*(subc)))
180 #define NBNXN_CI_DO_COUL(subc) (1<<(9+3*(subc)))
182 /* Cluster-pair Interaction masks
183 * Bit i*j-cluster-size + j tells if atom i and j interact.
185 // TODO: Rename according to convention when moving into Nbnxn namespace
186 /* All interaction mask is the same for all kernels */
187 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL = 0xffffffffU;
188 /* 4x4 kernel diagonal mask */
189 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG = 0x08ceU;
190 /* 4x2 kernel diagonal masks */
191 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
192 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
193 /* 4x8 kernel diagonal masks */
194 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
195 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
197 /* Simple pair-list i-unit */
200 int ci; /* i-cluster */
201 int shift; /* Shift vector index plus possible flags, see above */
202 int cj_ind_start; /* Start index into cj */
203 int cj_ind_end; /* End index into cj */
206 /* Grouped pair-list i-unit */
208 /* Returns the number of j-cluster groups in this entry */
209 int numJClusterGroups() const
211 return cj4_ind_end - cj4_ind_start;
214 int sci; /* i-super-cluster */
215 int shift; /* Shift vector index plus possible flags */
216 int cj4_ind_start; /* Start index into cj4 */
217 int cj4_ind_end; /* End index into cj4 */
220 /* Interaction data for a j-group for one warp */
223 // The i-cluster interactions mask for 1 warp
224 unsigned int imask = 0U;
225 // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
230 int cj[c_nbnxnGpuJgroupSize]; /* The 4 j-clusters */
231 nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data for 2 warps */
234 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
237 /* Constructor, sets no exclusions, so all atom pairs interacting */
240 for (unsigned int &pairEntry : pair)
242 pairEntry = NBNXN_INTERACTION_MASK_ALL;
246 /* Topology exclusion interaction bits per warp */
247 unsigned int pair[c_nbnxnGpuExclSize];
250 /* Cluster pairlist type for use on CPUs */
251 struct NbnxnPairlistCpu
253 gmx_cache_protect_t cp0;
255 int na_ci; /* The number of atoms per i-cluster */
256 int na_cj; /* The number of atoms per j-cluster */
257 real rlist; /* The radius for constructing the list */
258 FastVector<nbnxn_ci_t> ci; /* The i-cluster list */
259 FastVector<nbnxn_ci_t> ciOuter; /* The outer, unpruned i-cluster list */
261 FastVector<nbnxn_cj_t> cj; /* The j-cluster list, size ncj */
262 FastVector<nbnxn_cj_t> cjOuter; /* The outer, unpruned j-cluster list */
263 int ncjInUse; /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
265 int nci_tot; /* The total number of i clusters */
267 NbnxnPairlistCpuWork *work;
269 gmx_cache_protect_t cp1;
272 /* Cluster pairlist type, with extra hierarchies, for on the GPU
274 * NOTE: for better performance when combining lists over threads,
275 * all vectors should use default initialization. But when
276 * changing this, excl should be intialized when adding entries.
278 struct NbnxnPairlistGpu
282 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
284 NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
286 gmx_cache_protect_t cp0;
288 int na_ci; /* The number of atoms per i-cluster */
289 int na_cj; /* The number of atoms per j-cluster */
290 int na_sc; /* The number of atoms per super cluster */
291 real rlist; /* The radius for constructing the list */
292 // The i-super-cluster list, indexes into cj4;
293 gmx::HostVector<nbnxn_sci_t> sci;
294 // The list of 4*j-cluster groups
295 gmx::HostVector<nbnxn_cj4_t> cj4;
296 // Atom interaction bits (non-exclusions)
297 gmx::HostVector<nbnxn_excl_t> excl;
298 // The total number of i-clusters
301 NbnxnPairlistGpuWork *work;
303 gmx_cache_protect_t cp1;
306 struct nbnxn_pairlist_set_t
308 nbnxn_pairlist_set_t(const NbnxnListParameters &listParams);
310 int nnbl; /* number of lists */
311 NbnxnPairlistCpu **nbl; /* lists for CPU */
312 NbnxnPairlistCpu **nbl_work; /* work space for rebalancing lists */
313 NbnxnPairlistGpu **nblGpu; /* lists for GPU */
314 const NbnxnListParameters ¶ms; /* Pairlist parameters desribing setup and ranges */
315 gmx_bool bCombined; /* TRUE if lists get combined into one (the 1st) */
316 gmx_bool bSimple; /* TRUE if the list of of type "simple"
317 (na_sc=na_s, no super-clusters used) */
319 /* Counts for debug printing */
320 int natpair_ljq; /* Total number of atom pairs for LJ+Q kernel */
321 int natpair_lj; /* Total number of atom pairs for LJ kernel */
322 int natpair_q; /* Total number of atom pairs for Q kernel */
323 std::vector<t_nblist *> nbl_fep; /* List of free-energy atom pair interactions */
324 int64_t outerListCreationStep; /* Step at which the outer list was created */
328 nbatXYZ, nbatXYZQ, nbatX4, nbatX8
331 // Struct that holds force and energy output buffers
332 struct nbnxn_atomdata_output_t
336 * \param[in] kernelType Type of non-bonded kernel
337 * \param[in] numEnergyGroups The number of energy groups
338 * \param[in] simdEnergyBufferStride Stride for entries in the energy buffers for SIMD kernels
339 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
341 nbnxn_atomdata_output_t(Nbnxm::KernelType kernelType,
343 int simdEnergyBUfferStride,
344 gmx::PinningPolicy pinningPolicy);
346 gmx::HostVector<real> f; // f, size natoms*fstride
347 gmx::HostVector<real> fshift; // Shift force array, size SHIFTS*DIM
348 gmx::HostVector<real> Vvdw; // Temporary Van der Waals group energy storage
349 gmx::HostVector<real> Vc; // Temporary Coulomb group energy storage
350 AlignedVector<real> VSvdw; // Temporary SIMD Van der Waals group energy storage
351 AlignedVector<real> VSc; // Temporary SIMD Coulomb group energy storage
354 /* Block size in atoms for the non-bonded thread force-buffer reduction,
355 * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
356 * Should be small to reduce the reduction and zeroing cost,
357 * but too small will result in overhead.
358 * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
361 #define NBNXN_BUFFERFLAG_SIZE 8
363 #define NBNXN_BUFFERFLAG_SIZE 16
366 /* We store the reduction flags as gmx_bitmask_t.
367 * This limits the number of flags to BITMASK_SIZE.
369 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
371 /* Flags for telling if threads write to force output buffers */
373 int nflag; /* The number of flag blocks */
374 gmx_bitmask_t *flag; /* Bit i is set when thread i writes to a cell-block */
375 int flag_nalloc; /* Allocation size of cxy_flag */
376 } nbnxn_buffer_flags_t;
378 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
380 ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
383 /* Struct that stores atom related data for the nbnxn module
385 * Note: performance would improve slightly when all std::vector containers
386 * in this struct would not initialize during resize().
388 struct nbnxn_atomdata_t
389 { //NOLINT(clang-analyzer-optin.performance.Padding)
394 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
396 Params(gmx::PinningPolicy pinningPolicy);
398 // The number of different atom types
400 // Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
401 gmx::HostVector<real> nbfp;
402 // Combination rule, see enum defined above
404 // LJ parameters per atom type, size numTypes*2
405 gmx::HostVector<real> nbfp_comb;
406 // As nbfp, but with a stride for the present SIMD architecture
407 AlignedVector<real> nbfp_aligned;
408 // Atom types per atom
409 gmx::HostVector<int> type;
410 // LJ parameters per atom for fast SIMD loading
411 gmx::HostVector<real> lj_comb;
412 // Charges per atom, not set with format nbatXYZQ
413 gmx::HostVector<real> q;
414 // The number of energy groups
418 // The energy groups, one int entry per cluster, only set when needed
419 gmx::HostVector<int> energrp;
422 // Diagonal and topology exclusion helper data for all SIMD kernels
427 // Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
428 AlignedVector<real> diagonal_4xn_j_minus_i;
429 // Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
430 AlignedVector<real> diagonal_2xnn_j_minus_i;
431 // Filters for topology exclusion masks for the SIMD kernels
432 AlignedVector<uint32_t> exclusion_filter;
433 // Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
434 AlignedVector<uint64_t> exclusion_filter64;
435 // Array of masks needed for exclusions
436 AlignedVector<real> interaction_array;
441 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
443 nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
445 /* Returns a const reference to the parameters */
446 const Params ¶ms() const
451 /* Returns a non-const reference to the parameters */
452 Params ¶msDeprecated()
457 /* Returns the current total number of atoms stored */
463 /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
464 gmx::ArrayRef<const real> x() const
469 /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
470 gmx::ArrayRef<real> x()
475 /* Resizes the coordinate buffer and sets the number of atoms */
476 void resizeCoordinateBuffer(int numAtoms);
478 /* Resizes the force buffers for the current number of atoms */
479 void resizeForceBuffers();
482 // The LJ and charge parameters
484 // The total number of atoms currently stored
487 int natoms_local; /* Number of local atoms */
488 int XFormat; /* The format of x (and q), enum */
489 int FFormat; /* The format of f, enum */
490 gmx_bool bDynamicBox; /* Do we need to update shift_vec every step? */
491 gmx::HostVector<gmx::RVec> shift_vec; /* Shift vectors, copied from t_forcerec */
492 int xstride; /* stride for a coordinate in x (usually 3 or 4) */
493 int fstride; /* stride for a coordinate in f (usually 3 or 4) */
495 gmx::HostVector<real> x_; /* x and possibly q, size natoms*xstride */
498 // Masks for handling exclusions in the SIMD kernels
499 const SimdMasks simdMasks;
502 std::vector<nbnxn_atomdata_output_t> out; /* Output data structures, 1 per thread */
504 /* Reduction related data */
505 gmx_bool bUseBufferFlags; /* Use the flags or operate on all atoms */
506 nbnxn_buffer_flags_t buffer_flags; /* Flags for buffer zeroing+reduc. */
507 gmx_bool bUseTreeReduce; /* Use tree for force reduction */
508 tMPI_Atomic *syncStep; /* Synchronization step for tree reduce */