2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
37 * \brief OpenCL pruning kernel.
39 * OpenCL 1.2 support is expected; tested on AMD GCN and NVIDIA CC >3.0.
41 * \author Szilárd Páll <pall.szilard@gmail.com>
42 * \ingroup module_nbnxm
45 #include "nbnxm_ocl_kernel_utils.clh"
47 /* Note: the AMD compiler testing was done with (fglrx 15.12) performs best with wg
48 * size 256 (this is an artificial compiler limitation). The compiler is also
49 * sensitive to tidx/widx declaration and warp_any initialization.
50 * With the current tweaks the regular prune kenel achieves 90%, the rolling 100%
51 * occupancy with both Fiji and Hawaii.
52 * TODO: if the wg size limit is removed in an upcoming AMD compiler the NTHREAD_Z=4
53 * should be revisited.
56 #define NTHREAD_Z GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
58 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, NTHREAD_Z)))
59 #ifdef HAVE_FRESH_LIST
60 __kernel void nbnxn_kernel_prune_opencl
62 __kernel void nbnxn_kernel_prune_rolling_opencl
65 cl_nbparam_params_t nbparam_params,
66 const __global float4 *restrict xq,
67 const __global float *restrict shift_vec,
68 const __global nbnxn_sci_t *pl_sci,
69 __global nbnxn_cj4_t *pl_cj4,
70 #if !defined HAVE_FRESH_LIST
73 __global unsigned int *restrict prePrunedImask,
79 /* convenience variables */
80 const cl_nbparam_params_t *const nbparam = &nbparam_params;
82 const float rlistOuter_sq = nbparam->rlistOuter_sq;
83 const float rlistInner_sq = nbparam->rlistInner_sq;
85 /* thread/block/warp id-s */
86 const unsigned int tidxi = get_local_id(0);
87 const unsigned int tidxj = get_local_id(1);
88 const unsigned int tidx = get_local_id(1) * get_local_size(0) + get_local_id(0);
90 const unsigned int tidxz = 0;
92 const unsigned int tidxz = get_local_id(2);
94 const unsigned int bidx = get_group_id(0);
95 const unsigned int widx = tidx / WARP_SIZE;
97 #ifdef HAVE_FRESH_LIST
98 const bool haveFreshList = true;
100 const bool haveFreshList = false;
103 // TODO move these consts to utils and unify their use with the nonbonded kernels
104 const int c_numClPerSupercl = NCL_PER_SUPERCL;
105 const int c_clSize = CL_SIZE;
107 // TODO pass this value at compile-time as a macro
108 const int c_nbnxnGpuClusterpairSplit = 2;
110 /*! i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
111 const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
113 #define LOCAL_OFFSET (xib + c_numClPerSupercl * c_clSize)
114 /* shmem buffer for i cj pre-loading */
117 cjs = (((__local int *)(LOCAL_OFFSET)) + tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize);
119 /* Offset calculated using xib because cjs depends on on tidxz! */
120 #define LOCAL_OFFSET (((__local int *)(xib + c_numClPerSupercl * c_clSize)) + (NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize))
122 #if !USE_SUBGROUP_ANY
123 /* Local buffer used to implement __any warp vote function from CUDA.
124 volatile is used to avoid compiler optimizations for AMD builds. */
125 volatile __local uint *const warp_any = (__local uint*)(LOCAL_OFFSET);
126 const unsigned int warpVoteSlot = NTHREAD_Z*tidxz + widx;
127 /* Initialise warp vote.*/
128 if (tidx == 0 || tidx == 32)
130 warp_any[warpVoteSlot] = 0;
133 __local uint *const warp_any = 0;
134 const unsigned int warpVoteSlot = 0;
138 const nbnxn_sci_t nb_sci = pl_sci[bidx*numParts + part]; /* my i super-cluster's index = sciOffset + current bidx * numParts + part */
139 const int sci = nb_sci.sci; /* super-cluster */
140 const int cij4_start = nb_sci.cj4_ind_start; /* first ...*/
141 const int cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
145 for (int i = 0; i < NCL_PER_SUPERCL; i += CL_SIZE)
147 /* Pre-load i-atom x and q into shared memory */
148 const int ci = sci * c_numClPerSupercl + tidxj+i;
149 const int ai = ci * c_clSize + tidxi;
151 /* We don't need q, but using float4 in shmem avoids bank conflicts */
152 const float4 tmp = xq[ai];
153 const float4 xi = tmp + (float4)(shift_vec[3 * nb_sci.shift], shift_vec[3 * nb_sci.shift + 1], shift_vec[3 * nb_sci.shift + 2], 0.0f);
154 xib[(tidxj + i) * c_clSize + tidxi] = xi;
157 barrier(CLK_LOCAL_MEM_FENCE);
160 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
161 for (int j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
163 unsigned int imaskFull, imaskCheck, imaskNew;
167 /* Read the mask from the list transferred from the CPU */
168 imaskFull = pl_cj4[j4].imei[widx].imask;
169 /* We attempt to prune all pairs present in the original list */
170 imaskCheck = imaskFull;
175 /* Read the mask from the "warp-pruned" by rlistOuter mask array */
176 imaskFull = prePrunedImask[j4*c_nbnxnGpuClusterpairSplit + widx];
177 /* Read the old rolling pruned mask, use as a base for new */
178 imaskNew = pl_cj4[j4].imei[widx].imask;
179 /* We only need to check pairs with different mask */
180 imaskCheck = (imaskNew ^ imaskFull);
183 preloadCj4(&cjs, pl_cj4[j4].cj, tidxi, tidxj, imaskCheck != 0u);
188 for (int jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
190 if (imaskCheck & (superClInteractionMask << (jm * c_numClPerSupercl)))
192 unsigned int mask_ji = (1U << (jm * c_numClPerSupercl));
194 const int cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
195 const int aj = cj * c_clSize + tidxj;
197 /* load j atom data */
198 const float4 tmp = xq[aj];
199 const float3 xj = (float3)(tmp.xyz);
202 for (int i = 0; i < c_numClPerSupercl; i++)
204 if (imaskCheck & mask_ji)
206 /* load i-cluster coordinates from shmem */
207 const float4 xi = xib[i * c_clSize + tidxi];
209 /* distance between i and j atoms */
210 const float3 rv = (float3)(xi.xyz) - xj;
211 const float r2 = norm2(rv);
215 /* If _none_ of the atoms pairs are in cutoff range,
216 the bit corresponding to the current
217 cluster-pair in imask gets set to 0. */
218 if (!gmx_sub_group_any(warp_any, warpVoteSlot, r2 < rlistOuter_sq))
220 imaskFull &= ~mask_ji;
223 /* If any atom pair is within range, set the bit
224 corresponding to the current cluster-pair. */
225 if (gmx_sub_group_any(warp_any, warpVoteSlot, r2 < rlistInner_sq))
231 /* shift the mask bit by 1 */
237 #ifdef HAVE_FRESH_LIST
239 /* copy the list pruned to rlistOuter to a separate buffer */
240 prePrunedImask[j4*c_nbnxnGpuClusterpairSplit + widx] = imaskFull;
243 /* update the imask with only the pairs up to rlistInner */
244 pl_cj4[j4].imei[widx].imask = imaskNew;
249 #undef USE_CJ_PREFETCH