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/real.h"
51 // This file with constants is separate from this file to be able
52 // to include it during OpenCL jitting without including config.h
53 #include "gromacs/nbnxm/constants.h"
55 struct NbnxnPairlistCpuWork;
56 struct NbnxnPairlistGpuWork;
59 /* Convenience type for vector with aligned memory */
61 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
63 /* Convenience type for vector that avoids initialization at resize() */
65 using FastVector = std::vector < T, gmx::DefaultInitializationAllocator < T>>;
69 /*! \brief The setup for generating and pruning the nbnxn pair list.
71 * Without dynamic pruning rlistOuter=rlistInner.
73 struct NbnxnListParameters
75 /*! \brief Constructor producing a struct with dynamic pruning disabled
77 NbnxnListParameters(real rlist) :
78 useDynamicPruning(false),
86 bool useDynamicPruning; //!< Are we using dynamic pair-list pruning
87 int nstlistPrune; //!< Pair-list dynamic pruning interval
88 real rlistOuter; //!< Cut-off of the larger, outer pair-list
89 real rlistInner; //!< Cut-off of the smaller, inner pair-list
90 int numRollingParts; //!< The number parts to divide the pair-list into for rolling pruning, a value of 1 gives no rolling pruning
95 /* With CPU kernels the i-cluster size is always 4 atoms. */
96 static constexpr int c_nbnxnCpuIClusterSize = 4;
98 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
99 #if GMX_GPU == GMX_GPU_OPENCL
100 static constexpr int c_nbnxnGpuClusterSize = GMX_OPENCL_NB_CLUSTER_SIZE;
102 static constexpr int c_nbnxnGpuClusterSize = 8;
105 /* The number of clusters in a pair-search cell, used for GPU */
106 static constexpr int c_gpuNumClusterPerCellZ = 2;
107 static constexpr int c_gpuNumClusterPerCellY = 2;
108 static constexpr int c_gpuNumClusterPerCellX = 2;
109 static constexpr int c_gpuNumClusterPerCell = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
112 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
113 * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
115 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
117 /* The fixed size of the exclusion mask array for a half cluster pair */
118 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
120 /* A buffer data structure of 64 bytes
121 * to be placed at the beginning and end of structs
122 * to avoid cache invalidation of the real contents
123 * of the struct by writes to neighboring memory.
127 } gmx_cache_protect_t;
129 /* Abstract type for pair searching data */
130 typedef struct nbnxn_search * nbnxn_search_t;
132 /* Function that should return a pointer *ptr to memory
134 * Error handling should be done within this function.
136 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
138 /* Function that should free the memory pointed to by *ptr.
139 * NULL should not be passed to this function.
141 typedef void nbnxn_free_t (void *ptr);
143 /* This is the actual cluster-pair list j-entry.
144 * cj is the j-cluster.
145 * The interaction bits in excl are indexed i-major, j-minor.
146 * The cj entries are sorted such that ones with exclusions come first.
147 * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
148 * is found, all subsequent j-entries in the i-entry also have full masks.
152 int cj; /* The j-cluster */
153 unsigned int excl; /* The exclusion (interaction) bits */
156 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
157 * The upper bits contain information for non-bonded kernel optimization.
158 * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
159 * But three flags can be used to skip interactions, currently only for subc=0
160 * !(shift & NBNXN_CI_DO_LJ(subc)) => we can skip LJ for all pairs
161 * shift & NBNXN_CI_HALF_LJ(subc) => we can skip LJ for the second half of i
162 * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
164 #define NBNXN_CI_SHIFT 127
165 #define NBNXN_CI_DO_LJ(subc) (1<<(7+3*(subc)))
166 #define NBNXN_CI_HALF_LJ(subc) (1<<(8+3*(subc)))
167 #define NBNXN_CI_DO_COUL(subc) (1<<(9+3*(subc)))
169 /* Cluster-pair Interaction masks
170 * Bit i*j-cluster-size + j tells if atom i and j interact.
172 // TODO: Rename according to convention when moving into Nbnxn namespace
173 /* All interaction mask is the same for all kernels */
174 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL = 0xffffffffU;
175 /* 4x4 kernel diagonal mask */
176 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG = 0x08ceU;
177 /* 4x2 kernel diagonal masks */
178 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
179 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
180 /* 4x8 kernel diagonal masks */
181 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
182 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
184 /* Simple pair-list i-unit */
187 int ci; /* i-cluster */
188 int shift; /* Shift vector index plus possible flags, see above */
189 int cj_ind_start; /* Start index into cj */
190 int cj_ind_end; /* End index into cj */
193 /* Grouped pair-list i-unit */
195 /* Returns the number of j-cluster groups in this entry */
196 int numJClusterGroups() const
198 return cj4_ind_end - cj4_ind_start;
201 int sci; /* i-super-cluster */
202 int shift; /* Shift vector index plus possible flags */
203 int cj4_ind_start; /* Start index into cj4 */
204 int cj4_ind_end; /* End index into cj4 */
207 /* Interaction data for a j-group for one warp */
210 // The i-cluster interactions mask for 1 warp
211 unsigned int imask = 0U;
212 // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
217 int cj[c_nbnxnGpuJgroupSize]; /* The 4 j-clusters */
218 nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data for 2 warps */
221 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
224 /* Constructor, sets no exclusions, so all atom pairs interacting */
227 for (unsigned int &pairEntry : pair)
229 pairEntry = NBNXN_INTERACTION_MASK_ALL;
233 /* Topology exclusion interaction bits per warp */
234 unsigned int pair[c_nbnxnGpuExclSize];
237 /* Cluster pairlist type for use on CPUs */
238 struct NbnxnPairlistCpu
240 gmx_cache_protect_t cp0;
242 int na_ci; /* The number of atoms per i-cluster */
243 int na_cj; /* The number of atoms per j-cluster */
244 real rlist; /* The radius for constructing the list */
245 FastVector<nbnxn_ci_t> ci; /* The i-cluster list */
246 FastVector<nbnxn_ci_t> ciOuter; /* The outer, unpruned i-cluster list */
248 FastVector<nbnxn_cj_t> cj; /* The j-cluster list, size ncj */
249 FastVector<nbnxn_cj_t> cjOuter; /* The outer, unpruned j-cluster list */
250 int ncjInUse; /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
252 int nci_tot; /* The total number of i clusters */
254 NbnxnPairlistCpuWork *work;
256 gmx_cache_protect_t cp1;
259 /* Cluster pairlist type, with extra hierarchies, for on the GPU
261 * NOTE: for better performance when combining lists over threads,
262 * all vectors should use default initialization. But when
263 * changing this, excl should be intialized when adding entries.
265 struct NbnxnPairlistGpu
269 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
271 NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
273 gmx_cache_protect_t cp0;
275 int na_ci; /* The number of atoms per i-cluster */
276 int na_cj; /* The number of atoms per j-cluster */
277 int na_sc; /* The number of atoms per super cluster */
278 real rlist; /* The radius for constructing the list */
279 // The i-super-cluster list, indexes into cj4;
280 gmx::HostVector<nbnxn_sci_t> sci;
281 // The list of 4*j-cluster groups
282 gmx::HostVector<nbnxn_cj4_t> cj4;
283 // Atom interaction bits (non-exclusions)
284 gmx::HostVector<nbnxn_excl_t> excl;
285 // The total number of i-clusters
288 NbnxnPairlistGpuWork *work;
290 gmx_cache_protect_t cp1;
293 struct nbnxn_pairlist_set_t
295 int nnbl; /* number of lists */
296 NbnxnPairlistCpu **nbl; /* lists for CPU */
297 NbnxnPairlistCpu **nbl_work; /* work space for rebalancing lists */
298 NbnxnPairlistGpu **nblGpu; /* lists for GPU */
299 gmx_bool bCombined; /* TRUE if lists get combined into one (the 1st) */
300 gmx_bool bSimple; /* TRUE if the list of of type "simple"
301 (na_sc=na_s, no super-clusters used) */
302 int natpair_ljq; /* Total number of atom pairs for LJ+Q kernel */
303 int natpair_lj; /* Total number of atom pairs for LJ kernel */
304 int natpair_q; /* Total number of atom pairs for Q kernel */
305 t_nblist **nbl_fep; /* List of free-energy atom pair interactions */
306 int64_t outerListCreationStep; /* Step at which the outer list was created */
310 nbatXYZ, nbatXYZQ, nbatX4, nbatX8
313 // Struct that holds force and energy output buffers
314 struct nbnxn_atomdata_output_t
318 * \param[in] nb_kernel_type Type of non-bonded kernel
319 * \param[in] numEnergyGroups The number of energy groups
320 * \param[in] simdEnergyBufferStride Stride for entries in the energy buffers for SIMD kernels
321 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
323 nbnxn_atomdata_output_t(int nb_kernel_type,
325 int simdEnergyBUfferStride,
326 gmx::PinningPolicy pinningPolicy);
328 gmx::HostVector<real> f; // f, size natoms*fstride
329 gmx::HostVector<real> fshift; // Shift force array, size SHIFTS*DIM
330 gmx::HostVector<real> Vvdw; // Temporary Van der Waals group energy storage
331 gmx::HostVector<real> Vc; // Temporary Coulomb group energy storage
332 AlignedVector<real> VSvdw; // Temporary SIMD Van der Waals group energy storage
333 AlignedVector<real> VSc; // Temporary SIMD Coulomb group energy storage
336 /* Block size in atoms for the non-bonded thread force-buffer reduction,
337 * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
338 * Should be small to reduce the reduction and zeroing cost,
339 * but too small will result in overhead.
340 * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
343 #define NBNXN_BUFFERFLAG_SIZE 8
345 #define NBNXN_BUFFERFLAG_SIZE 16
348 /* We store the reduction flags as gmx_bitmask_t.
349 * This limits the number of flags to BITMASK_SIZE.
351 #define NBNXN_BUFFERFLAG_MAX_THREADS (BITMASK_SIZE)
353 /* Flags for telling if threads write to force output buffers */
355 int nflag; /* The number of flag blocks */
356 gmx_bitmask_t *flag; /* Bit i is set when thread i writes to a cell-block */
357 int flag_nalloc; /* Allocation size of cxy_flag */
358 } nbnxn_buffer_flags_t;
360 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
362 ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
365 /* Struct that stores atom related data for the nbnxn module
367 * Note: performance would improve slightly when all std::vector containers
368 * in this struct would not initialize during resize().
370 struct nbnxn_atomdata_t
371 { //NOLINT(clang-analyzer-optin.performance.Padding)
376 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
378 Params(gmx::PinningPolicy pinningPolicy);
380 // The number of different atom types
382 // Lennard-Jone 6*C6 and 12*C12 parameters, size numTypes*2*2
383 gmx::HostVector<real> nbfp;
384 // Combination rule, see enum defined above
386 // LJ parameters per atom type, size numTypes*2
387 gmx::HostVector<real> nbfp_comb;
388 // As nbfp, but with a stride for the present SIMD architecture
389 AlignedVector<real> nbfp_aligned;
390 // Atom types per atom
391 gmx::HostVector<int> type;
392 // LJ parameters per atom for fast SIMD loading
393 gmx::HostVector<real> lj_comb;
394 // Charges per atom, not set with format nbatXYZQ
395 gmx::HostVector<real> q;
396 // The number of energy groups
400 // The energy groups, one int entry per cluster, only set when needed
401 gmx::HostVector<int> energrp;
404 // Diagonal and topology exclusion helper data for all SIMD kernels
409 // Helper data for setting up diagonal exclusion masks in the SIMD 4xN kernels
410 AlignedVector<real> diagonal_4xn_j_minus_i;
411 // Helper data for setting up diaginal exclusion masks in the SIMD 2xNN kernels
412 AlignedVector<real> diagonal_2xnn_j_minus_i;
413 // Filters for topology exclusion masks for the SIMD kernels
414 AlignedVector<uint32_t> exclusion_filter;
415 // Filters for topology exclusion masks for double SIMD kernels without SIMD int32 logical support
416 AlignedVector<uint64_t> exclusion_filter64;
417 // Array of masks needed for exclusions
418 AlignedVector<real> interaction_array;
423 * \param[in] pinningPolicy Sets the pinning policy for all data that might be transfered to a GPU
425 nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy);
427 /* Returns a const reference to the parameters */
428 const Params ¶ms() const
433 /* Returns a non-const reference to the parameters */
434 Params ¶msDeprecated()
439 /* Returns the current total number of atoms stored */
445 /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
446 gmx::ArrayRef<const real> x() const
451 /* Return the coordinate buffer, and q with xFormat==nbatXYZQ */
452 gmx::ArrayRef<real> x()
457 /* Resizes the coordinate buffer and sets the number of atoms */
458 void resizeCoordinateBuffer(int numAtoms);
460 /* Resizes the force buffers for the current number of atoms */
461 void resizeForceBuffers();
464 // The LJ and charge parameters
466 // The total number of atoms currently stored
469 int natoms_local; /* Number of local atoms */
470 int XFormat; /* The format of x (and q), enum */
471 int FFormat; /* The format of f, enum */
472 gmx_bool bDynamicBox; /* Do we need to update shift_vec every step? */
473 gmx::HostVector<gmx::RVec> shift_vec; /* Shift vectors, copied from t_forcerec */
474 int xstride; /* stride for a coordinate in x (usually 3 or 4) */
475 int fstride; /* stride for a coordinate in f (usually 3 or 4) */
477 gmx::HostVector<real> x_; /* x and possibly q, size natoms*xstride */
480 // Masks for handling exclusions in the SIMD kernels
481 const SimdMasks simdMasks;
484 std::vector<nbnxn_atomdata_output_t> out; /* Output data structures, 1 per thread */
486 /* Reduction related data */
487 gmx_bool bUseBufferFlags; /* Use the flags or operate on all atoms */
488 nbnxn_buffer_flags_t buffer_flags; /* Flags for buffer zeroing+reduc. */
489 gmx_bool bUseTreeReduce; /* Use tree for force reduction */
490 tMPI_Atomic *syncStep; /* Synchronization step for tree reduce */