Report name of missing OpenCL kernels
[gromacs.git] / src / gromacs / nbnxm / opencl / nbnxm_ocl_kernel_pruneonly.clh
blob16dd0c4962dbf0b853403b9ec955b6029b1d4e16
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
36 /*! \internal \file
37  *  \brief OpenCL pruning kernel.
38  *
39  *  OpenCL 1.2 support is expected; tested on AMD GCN and NVIDIA CC >3.0.
40  *
41  *  \author Szilárd Páll <pall.szilard@gmail.com>
42  *  \ingroup module_nbnxm
43  */
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.
54  *
55  */
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
61 nbnxn_kernel_prune_opencl
62 #else
63 __kernel void
64 nbnxn_kernel_prune_rolling_opencl
65 #endif
66         (cl_nbparam_params_t nbparam_params,
67          const __global float4* restrict xq,
68          const __global float* restrict shift_vec,
69          const __global nbnxn_sci_t* pl_sci,
70          __global nbnxn_cj4_t* pl_cj4,
71 #if !defined HAVE_FRESH_LIST
72          const
73 #endif
74          __global unsigned int* restrict prePrunedImask,
75          int                             numParts,
76          int                             part,
77          __local float4* xib)
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 int tidxi = get_local_id(0);
87     const int tidxj = get_local_id(1);
88     const int tidx  = (int)(get_local_id(1) * get_local_size(0) + get_local_id(0));
89 #if NTHREAD_Z == 1
90     const int tidxz = 0;
91 #else
92     const int          tidxz         = get_local_id(2);
93 #endif
94     const int bidx = get_group_id(0);
95     const int widx = tidx / WARP_SIZE;
97 #ifdef HAVE_FRESH_LIST
98     const bool haveFreshList = true;
99 #else
100     const bool         haveFreshList = false;
101 #endif
103     // TODO move these consts to utils and unify their use with the nonbonded kernels
104     const int c_clSize = CL_SIZE;
106     // TODO pass this value at compile-time as a macro
107     const int c_nbnxnGpuClusterpairSplit = 2;
109     /*! i-cluster interaction mask for a super-cluster with all c_nbnxnGpuNumClusterPerSupercluster=8 bits set */
110     const unsigned superClInteractionMask = ((1U << c_nbnxnGpuNumClusterPerSupercluster) - 1U);
112 #define LOCAL_OFFSET (xib + c_nbnxnGpuNumClusterPerSupercluster * c_clSize)
113     /* shmem buffer for i cj pre-loading */
114     CjType cjs = 0;
115 #if USE_CJ_PREFETCH
116     cjs = (((__local int*)(LOCAL_OFFSET)) + tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize);
117 #    undef LOCAL_OFFSET
118 /* Offset calculated using xib because cjs depends on on tidxz! */
119 #    define LOCAL_OFFSET                                                        \
120         (((__local int*)(xib + c_nbnxnGpuNumClusterPerSupercluster * c_clSize)) \
121          + (NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize))
122 #endif
123 #if !USE_SUBGROUP_ANY
124     /* Local buffer used to implement __any warp vote function from CUDA.
125        volatile is used to avoid compiler optimizations for AMD builds. */
126     //NOLINTNEXTLINE(google-readability-casting)
127     volatile __local int* const warp_any     = (__local int*)(LOCAL_OFFSET);
128     const int                   warpVoteSlot = NTHREAD_Z * tidxz + widx;
129     /* Initialise warp vote.*/
130     // if (tidx == 0 || tidx == 32)
131     if (tidx % (c_clSize * c_clSize / c_nbnxnGpuClusterpairSplit) == 0)
132     {
133         warp_any[warpVoteSlot] = 0;
134     }
135 #else
136     __local int* const warp_any      = 0;
137     const int          warpVoteSlot  = 0;
138 #endif
139 #undef LOCAL_OFFSET
141     const nbnxn_sci_t nb_sci =
142             pl_sci[bidx * numParts + part]; /* my i super-cluster's index = sciOffset + current bidx * numParts + part */
143     const int sci        = nb_sci.sci;           /* super-cluster */
144     const int cij4_start = nb_sci.cj4_ind_start; /* first ...*/
145     const int cij4_end   = nb_sci.cj4_ind_end;   /* and last index of j clusters */
147     if (tidxz == 0)
148     {
149         for (int i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i += CL_SIZE)
150         {
151             /* Pre-load i-atom x and q into shared memory */
152             const int ci = sci * c_nbnxnGpuNumClusterPerSupercluster + tidxj + i;
153             const int ai = ci * c_clSize + tidxi;
155             /* We don't need q, but using float4 in shmem avoids bank conflicts */
156             const float4 tmp = xq[ai];
157             const float4 xi  = tmp
158                               + (float4)(shift_vec[3 * nb_sci.shift], shift_vec[3 * nb_sci.shift + 1],
159                                          shift_vec[3 * nb_sci.shift + 2], 0.0F);
160             xib[(tidxj + i) * c_clSize + tidxi] = xi;
161         }
162     }
163     barrier(CLK_LOCAL_MEM_FENCE);
166     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
167     for (int j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
168     {
169         unsigned int imaskFull, imaskCheck, imaskNew;
171         if (haveFreshList)
172         {
173             /* Read the mask from the list transferred from the CPU */
174             imaskFull = pl_cj4[j4].imei[widx].imask;
175             /* We attempt to prune all pairs present in the original list */
176             imaskCheck = imaskFull;
177             imaskNew   = 0;
178         }
179         else
180         {
181             /* Read the mask from the "warp-pruned" by rlistOuter mask array */
182             imaskFull = prePrunedImask[j4 * c_nbnxnGpuClusterpairSplit + widx];
183             /* Read the old rolling pruned mask, use as a base for new */
184             imaskNew = pl_cj4[j4].imei[widx].imask;
185             /* We only need to check pairs with different mask */
186             imaskCheck = (imaskNew ^ imaskFull);
187         }
189         preloadCj4(&cjs, pl_cj4[j4].cj, tidxi, tidxj, imaskCheck != 0U);
191         if (imaskCheck)
192         {
193 #pragma unroll 4
194             for (int jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
195             {
196                 if (imaskCheck & (superClInteractionMask << (jm * c_nbnxnGpuNumClusterPerSupercluster)))
197                 {
198                     unsigned int mask_ji = (1U << (jm * c_nbnxnGpuNumClusterPerSupercluster));
200                     const int cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
201                     const int aj = cj * c_clSize + tidxj;
203                     /* load j atom data */
204                     const float4 tmp = xq[aj];
205                     const float3 xj  = (float3)(tmp.xyz);
207 #pragma unroll 8
208                     for (int i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
209                     {
210                         if (imaskCheck & mask_ji)
211                         {
212                             /* load i-cluster coordinates from shmem */
213                             const float4 xi = xib[i * c_clSize + tidxi];
215                             /* distance between i and j atoms */
216                             const float3 rv = (float3)(xi.xyz) - xj;
217                             const float  r2 = norm2(rv);
219                             if (haveFreshList)
220                             {
221                                 /* If _none_ of the atoms pairs are in cutoff range,
222                                    the bit corresponding to the current
223                                    cluster-pair in imask gets set to 0. */
224                                 if (!gmx_sub_group_any(warp_any, warpVoteSlot, r2 < rlistOuter_sq))
225                                 {
226                                     imaskFull &= ~mask_ji;
227                                 }
228                             }
229                             /* If any atom pair is within range, set the bit
230                                corresponding to the current cluster-pair. */
231                             if (gmx_sub_group_any(warp_any, warpVoteSlot, r2 < rlistInner_sq))
232                             {
233                                 imaskNew |= mask_ji;
234                             }
235                         }
237                         /* shift the mask bit by 1 */
238                         mask_ji += mask_ji;
239                     }
240                 }
241             }
243 #ifdef HAVE_FRESH_LIST
244             {
245                 /* copy the list pruned to rlistOuter to a separate buffer */
246                 prePrunedImask[j4 * c_nbnxnGpuClusterpairSplit + widx] = imaskFull;
247             }
248 #endif
249             /* update the imask with only the pairs up to rlistInner */
250             pl_cj4[j4].imei[widx].imask = imaskNew;
251         }
252     }
255 #undef USE_CJ_PREFETCH