2c1e9c196f21c45f6440455e10ce6b2e7abe3762
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_ocl / nbnxn_ocl_kernel_pruneonly.clh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018, 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  */
35
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_mdlib
43  */
44
45 #include "nbnxn_ocl_kernel_utils.clh"
46
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
57
58 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, NTHREAD_Z)))
59 #ifdef HAVE_FRESH_LIST
60 __kernel void nbnxn_kernel_prune_opencl
61 #else
62 __kernel void nbnxn_kernel_prune_rolling_opencl
63 #endif
64 (
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
71         const
72 #endif
73         __global unsigned int *restrict  prePrunedImask,
74         int                              numParts,
75         int                              part,
76         __local  float4                 *xib
77 )
78 {
79     /* convenience variables */
80     const cl_nbparam_params_t *const nbparam = &nbparam_params;
81
82     const float                      rlistOuter_sq = nbparam->rlistOuter_sq;
83     const float                      rlistInner_sq = nbparam->rlistInner_sq;
84
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);
89 #if NTHREAD_Z == 1
90     const unsigned int tidxz  = 0;
91 #else
92     const unsigned int tidxz  = get_local_id(2);
93 #endif
94     const unsigned int bidx   = get_group_id(0);
95     const unsigned int widx   = tidx / WARP_SIZE;
96
97 #ifdef HAVE_FRESH_LIST
98     const bool haveFreshList = true;
99 #else
100     const bool haveFreshList = false;
101 #endif
102
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;
106     const int c_nbnxnGpuJgroupSize = NBNXN_GPU_JGROUP_SIZE;
107
108     // TODO pass this value at compile-time as a macro
109     const int c_nbnxnGpuClusterpairSplit = 2;
110     const int c_splitClSize              = c_clSize/c_nbnxnGpuClusterpairSplit;
111
112     /*! i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
113     const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
114
115     #define LOCAL_OFFSET (xib + c_numClPerSupercl * c_clSize)
116     /* shmem buffer for i cj pre-loading */
117     CjType cjs;
118 #if USE_CJ_PREFETCH
119     cjs = (((__local int *)(LOCAL_OFFSET)) + tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize);
120     #undef LOCAL_OFFSET
121     /* Offset calculated using xib because cjs depends on on tidxz! */
122     #define LOCAL_OFFSET (((__local int *)(xib + c_numClPerSupercl * c_clSize)) + (NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize))
123 #endif
124 #if !USE_SUBGROUP_ANY
125     /* Local buffer used to implement __any warp vote function from CUDA.
126        volatile is used to avoid compiler optimizations for AMD builds. */
127     volatile __local uint *const warp_any     = (__local uint*)(LOCAL_OFFSET);
128     const unsigned int           warpVoteSlot = NTHREAD_Z*tidxz + widx;
129     /* Initialise warp vote.*/
130     if (tidx == 0 || tidx == 32)
131     {
132         warp_any[warpVoteSlot] = 0;
133     }
134 #else
135     __local uint *const warp_any     = 0;
136     const unsigned int  warpVoteSlot = 0;
137 #endif
138     #undef LOCAL_OFFSET
139
140     const nbnxn_sci_t nb_sci      = pl_sci[bidx*numParts + part]; /* my i super-cluster's index = sciOffset + current bidx * numParts + part */
141     const int         sci         = nb_sci.sci;                   /* super-cluster */
142     const int         cij4_start  = nb_sci.cj4_ind_start;         /* first ...*/
143     const int         cij4_end    = nb_sci.cj4_ind_end;           /* and last index of j clusters */
144
145     if (tidxz == 0)
146     {
147         for (int i = 0; i < NCL_PER_SUPERCL; i += CL_SIZE)
148         {
149             /* Pre-load i-atom x and q into shared memory */
150             const int ci = sci * c_numClPerSupercl + tidxj+i;
151             const int ai = ci * c_clSize + tidxi;
152
153             /* We don't need q, but using float4 in shmem avoids bank conflicts */
154             const float4 tmp = xq[ai];
155             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);
156             xib[(tidxj + i) * c_clSize + tidxi] = xi;
157         }
158     }
159     barrier(CLK_LOCAL_MEM_FENCE);
160
161
162     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
163     for (int j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
164     {
165         unsigned int imaskFull, imaskCheck, imaskNew;
166
167         if (haveFreshList)
168         {
169             /* Read the mask from the list transferred from the CPU */
170             imaskFull  = pl_cj4[j4].imei[widx].imask;
171             /* We attempt to prune all pairs present in the original list */
172             imaskCheck = imaskFull;
173             imaskNew   = 0;
174         }
175         else
176         {
177             /* Read the mask from the "warp-pruned" by rlistOuter mask array */
178             imaskFull  = prePrunedImask[j4*c_nbnxnGpuClusterpairSplit + widx];
179             /* Read the old rolling pruned mask, use as a base for new */
180             imaskNew   = pl_cj4[j4].imei[widx].imask;
181             /* We only need to check pairs with different mask */
182             imaskCheck = (imaskNew ^ imaskFull);
183         }
184
185         preloadCj4(&cjs, pl_cj4[j4].cj, tidxi, tidxj, imaskCheck != 0u);
186
187         if (imaskCheck)
188         {
189 #pragma unroll 4
190             for (int jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
191             {
192                 if (imaskCheck & (superClInteractionMask << (jm * c_numClPerSupercl)))
193                 {
194                     unsigned int       mask_ji = (1U << (jm * c_numClPerSupercl));
195
196                     const int          cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
197                     const int          aj = cj * c_clSize + tidxj;
198
199                     /* load j atom data */
200                     const float4 tmp  = xq[aj];
201                     const float3 xj   = (float3)(tmp.xyz);
202
203 #pragma unroll 8
204                     for (int i = 0; i < c_numClPerSupercl; i++)
205                     {
206                         if (imaskCheck & mask_ji)
207                         {
208                             /* load i-cluster coordinates from shmem */
209                             const float4 xi = xib[i * c_clSize + tidxi];
210
211                             /* distance between i and j atoms */
212                             const float3 rv = (float3)(xi.xyz) - xj;
213                             const float  r2 = norm2(rv);
214
215                             if (haveFreshList)
216                             {
217                                 /* If _none_ of the atoms pairs are in cutoff range,
218                                    the bit corresponding to the current
219                                    cluster-pair in imask gets set to 0. */
220                                 if (!gmx_sub_group_any(warp_any, warpVoteSlot, r2 < rlistOuter_sq))
221                                 {
222                                     imaskFull &= ~mask_ji;
223                                 }
224                             }
225                             /* If any atom pair is within range, set the bit
226                                corresponding to the current cluster-pair. */
227                             if (gmx_sub_group_any(warp_any, warpVoteSlot, r2 < rlistInner_sq))
228                             {
229                                 imaskNew |= mask_ji;
230                             }
231                         }
232
233                         /* shift the mask bit by 1 */
234                         mask_ji += mask_ji;
235                     }
236                 }
237             }
238
239 #ifdef HAVE_FRESH_LIST
240             {
241                 /* copy the list pruned to rlistOuter to a separate buffer */
242                 prePrunedImask[j4*c_nbnxnGpuClusterpairSplit + widx] = imaskFull;
243             }
244 #endif
245             /* update the imask with only the pairs up to rlistInner */
246             pl_cj4[j4].imei[widx].imask = imaskNew;
247         }
248     }
249 }
250
251 #undef USE_CJ_PREFETCH