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/defaultinitializationallocator.h"
48 #include "gromacs/utility/enumerationhelpers.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"
57 struct NbnxnPairlistCpuWork;
58 struct NbnxnPairlistGpuWork;
61 /* Convenience type for vector with aligned memory */
63 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
65 /* Convenience type for vector that avoids initialization at resize() */
67 using FastVector = std::vector < T, gmx::DefaultInitializationAllocator < T>>;
69 /* With CPU kernels the i-cluster size is always 4 atoms. */
70 static constexpr int c_nbnxnCpuIClusterSize = 4;
72 /* With GPU kernels the i and j cluster size is 8 atoms for CUDA and can be set at compile time for OpenCL */
73 #if GMX_GPU == GMX_GPU_OPENCL
74 static constexpr int c_nbnxnGpuClusterSize = GMX_OPENCL_NB_CLUSTER_SIZE;
76 static constexpr int c_nbnxnGpuClusterSize = 8;
79 /* The number of clusters in a pair-search cell, used for GPU */
80 static constexpr int c_gpuNumClusterPerCellZ = 2;
81 static constexpr int c_gpuNumClusterPerCellY = 2;
82 static constexpr int c_gpuNumClusterPerCellX = 2;
83 static constexpr int c_gpuNumClusterPerCell = c_gpuNumClusterPerCellZ*c_gpuNumClusterPerCellY*c_gpuNumClusterPerCellX;
86 /* In CUDA the number of threads in a warp is 32 and we have cluster pairs
87 * of 8*8=64 atoms, so it's convenient to store data for cluster pair halves.
89 static constexpr int c_nbnxnGpuClusterpairSplit = 2;
91 /* The fixed size of the exclusion mask array for a half cluster pair */
92 static constexpr int c_nbnxnGpuExclSize = c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
94 //! The available pair list types
95 enum class PairlistType : int
104 //! Gives the i-cluster size for each pairlist type
105 static constexpr gmx::EnumerationArray<PairlistType, int> IClusterSizePerListType =
107 c_nbnxnCpuIClusterSize,
108 c_nbnxnCpuIClusterSize,
109 c_nbnxnCpuIClusterSize,
110 c_nbnxnGpuClusterSize
112 //! Gives the j-cluster size for each pairlist type
113 static constexpr gmx::EnumerationArray<PairlistType, int> JClusterSizePerListType =
118 c_nbnxnGpuClusterSize
121 /* A buffer data structure of 64 bytes
122 * to be placed at the beginning and end of structs
123 * to avoid cache invalidation of the real contents
124 * of the struct by writes to neighboring memory.
128 } gmx_cache_protect_t;
130 /* This is the actual cluster-pair list j-entry.
131 * cj is the j-cluster.
132 * The interaction bits in excl are indexed i-major, j-minor.
133 * The cj entries are sorted such that ones with exclusions come first.
134 * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
135 * is found, all subsequent j-entries in the i-entry also have full masks.
139 int cj; /* The j-cluster */
140 unsigned int excl; /* The exclusion (interaction) bits */
143 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
144 * The upper bits contain information for non-bonded kernel optimization.
145 * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
146 * But three flags can be used to skip interactions, currently only for subc=0
147 * !(shift & NBNXN_CI_DO_LJ(subc)) => we can skip LJ for all pairs
148 * shift & NBNXN_CI_HALF_LJ(subc) => we can skip LJ for the second half of i
149 * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
151 #define NBNXN_CI_SHIFT 127
152 #define NBNXN_CI_DO_LJ(subc) (1<<(7+3*(subc)))
153 #define NBNXN_CI_HALF_LJ(subc) (1<<(8+3*(subc)))
154 #define NBNXN_CI_DO_COUL(subc) (1<<(9+3*(subc)))
156 /* Cluster-pair Interaction masks
157 * Bit i*j-cluster-size + j tells if atom i and j interact.
159 // TODO: Rename according to convention when moving into Nbnxn namespace
160 /* All interaction mask is the same for all kernels */
161 constexpr unsigned int NBNXN_INTERACTION_MASK_ALL = 0xffffffffU;
162 /* 4x4 kernel diagonal mask */
163 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG = 0x08ceU;
164 /* 4x2 kernel diagonal masks */
165 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_0 = 0x0002U;
166 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J2_1 = 0x002fU;
167 /* 4x8 kernel diagonal masks */
168 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_0 = 0xf0f8fcfeU;
169 constexpr unsigned int NBNXN_INTERACTION_MASK_DIAG_J8_1 = 0x0080c0e0U;
171 /* Simple pair-list i-unit */
174 int ci; /* i-cluster */
175 int shift; /* Shift vector index plus possible flags, see above */
176 int cj_ind_start; /* Start index into cj */
177 int cj_ind_end; /* End index into cj */
180 /* Grouped pair-list i-unit */
182 /* Returns the number of j-cluster groups in this entry */
183 int numJClusterGroups() const
185 return cj4_ind_end - cj4_ind_start;
188 int sci; /* i-super-cluster */
189 int shift; /* Shift vector index plus possible flags */
190 int cj4_ind_start; /* Start index into cj4 */
191 int cj4_ind_end; /* End index into cj4 */
194 /* Interaction data for a j-group for one warp */
197 // The i-cluster interactions mask for 1 warp
198 unsigned int imask = 0U;
199 // Index into the exclusion array for 1 warp, default index 0 which means no exclusions
204 int cj[c_nbnxnGpuJgroupSize]; /* The 4 j-clusters */
205 nbnxn_im_ei_t imei[c_nbnxnGpuClusterpairSplit]; /* The i-cluster mask data for 2 warps */
208 /* Struct for storing the atom-pair interaction bits for a cluster pair in a GPU pairlist */
211 /* Constructor, sets no exclusions, so all atom pairs interacting */
214 for (unsigned int &pairEntry : pair)
216 pairEntry = NBNXN_INTERACTION_MASK_ALL;
220 /* Topology exclusion interaction bits per warp */
221 unsigned int pair[c_nbnxnGpuExclSize];
224 /* Cluster pairlist type for use on CPUs */
225 struct NbnxnPairlistCpu
229 gmx_cache_protect_t cp0;
231 int na_ci; /* The number of atoms per i-cluster */
232 int na_cj; /* The number of atoms per j-cluster */
233 real rlist; /* The radius for constructing the list */
234 FastVector<nbnxn_ci_t> ci; /* The i-cluster list */
235 FastVector<nbnxn_ci_t> ciOuter; /* The outer, unpruned i-cluster list */
237 FastVector<nbnxn_cj_t> cj; /* The j-cluster list, size ncj */
238 FastVector<nbnxn_cj_t> cjOuter; /* The outer, unpruned j-cluster list */
239 int ncjInUse; /* The number of j-clusters that are used by ci entries in this list, will be <= cj.size() */
241 int nci_tot; /* The total number of i clusters */
243 /* Working data storage for list construction */
244 std::unique_ptr<NbnxnPairlistCpuWork> work;
246 gmx_cache_protect_t cp1;
249 /* Cluster pairlist type, with extra hierarchies, for on the GPU
251 * NOTE: for better performance when combining lists over threads,
252 * all vectors should use default initialization. But when
253 * changing this, excl should be intialized when adding entries.
255 struct NbnxnPairlistGpu
259 * \param[in] pinningPolicy Sets the pinning policy for all buffers used on the GPU
261 NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy);
263 gmx_cache_protect_t cp0;
265 int na_ci; /* The number of atoms per i-cluster */
266 int na_cj; /* The number of atoms per j-cluster */
267 int na_sc; /* The number of atoms per super cluster */
268 real rlist; /* The radius for constructing the list */
269 // The i-super-cluster list, indexes into cj4;
270 gmx::HostVector<nbnxn_sci_t> sci;
271 // The list of 4*j-cluster groups
272 gmx::HostVector<nbnxn_cj4_t> cj4;
273 // Atom interaction bits (non-exclusions)
274 gmx::HostVector<nbnxn_excl_t> excl;
275 // The total number of i-clusters
278 /* Working data storage for list construction */
279 std::unique_ptr<NbnxnPairlistGpuWork> work;
281 gmx_cache_protect_t cp1;
284 //! Initializes a free-energy pair-list
285 void nbnxn_init_pairlist_fep(t_nblist *nl);