1b7ed8478abe06464a0c4da06403cde0646d46fc
[alexxy/gromacs.git] / src / gromacs / nbnxm / opencl / nbnxm_ocl_kernel_utils.clh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2016,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  *  \brief
38  *  Utility constant and function declaration for the OpenCL non-bonded kernels.
39  *  This header should be included once at the top level, just before the
40  *  kernels are included (has to be preceded by nbnxn_ocl_types.h).
41  *
42  *  \author Szilárd Páll <pall.szilard@gmail.com>
43  *  \ingroup module_nbnxm
44  */
45
46 #define GMX_DOUBLE 0
47
48 #include "gromacs/gpu_utils/device_utils.clh"
49 #include "gromacs/gpu_utils/vectype_ops.clh"
50 #include "gromacs/pbcutil/ishift.h"
51
52 #include "nbnxm_ocl_consts.h"
53
54 #define CL_SIZE (c_nbnxnGpuClusterSize)
55
56 #define WARP_SIZE (CL_SIZE * CL_SIZE / 2) // Currently only c_nbnxnGpuClusterpairSplit=2 supported
57
58 #if defined _NVIDIA_SOURCE_ || defined _AMD_SOURCE_
59 /* Currently we enable CJ prefetch for AMD/NVIDIA and disable it for other vendors
60  * Note that this should precede the kernel_utils include.
61  */
62 #    define USE_CJ_PREFETCH 1
63 #else
64 #    define USE_CJ_PREFETCH 0
65 #endif
66
67 #if defined cl_intel_subgroups || defined cl_khr_subgroups \
68         || (defined __OPENCL_VERSION__ && __OPENCL_VERSION__ >= 210)
69 #    define HAVE_SUBGROUP 1
70 #else
71 #    define HAVE_SUBGROUP 0
72 #endif
73
74 #ifdef cl_intel_subgroups
75 #    define HAVE_INTEL_SUBGROUP 1
76 #else
77 #    define HAVE_INTEL_SUBGROUP 0
78 #endif
79
80 #if defined _INTEL_SOURCE_
81 #    define SUBGROUP_SIZE 8
82 #elif defined _AMD_SOURCE_
83 #    define SUBGROUP_SIZE 64
84 #else
85 #    define SUBGROUP_SIZE 32
86 #endif
87
88 #define REDUCE_SHUFFLE (HAVE_INTEL_SUBGROUP && CL_SIZE == 4 && SUBGROUP_SIZE == WARP_SIZE)
89 #define USE_SUBGROUP_ANY (HAVE_SUBGROUP && SUBGROUP_SIZE == WARP_SIZE)
90 #define USE_SUBGROUP_PRELOAD HAVE_INTEL_SUBGROUP
91
92 /* 1.0 / sqrt(M_PI) */
93 #define M_FLOAT_1_SQRTPI 0.564189583547756F
94
95 //-------------------
96
97 #ifndef NBNXN_OPENCL_KERNEL_UTILS_CLH
98 #    define NBNXN_OPENCL_KERNEL_UTILS_CLH
99
100 #    if CL_SIZE == 8
101 #        define WARP_SIZE_LOG2 (5)
102 #        define CL_SIZE_LOG2 (3)
103 #    elif CL_SIZE == 4
104 #        define WARP_SIZE_LOG2 (3)
105 #        define CL_SIZE_LOG2 (2)
106 #    else
107 #        error unsupported CL_SIZE
108 #    endif
109
110 #    define CL_SIZE_SQ (CL_SIZE * CL_SIZE)
111 #    define FBUF_STRIDE (CL_SIZE_SQ)
112
113 #    define ONE_SIXTH_F 0.16666667F
114 #    define ONE_TWELVETH_F 0.08333333F
115 #    define HALF_F 0.5F
116
117
118 #    ifdef __GNUC__
119 /* GCC, clang, and some ICC pretending to be GCC */
120 #        define gmx_unused __attribute__((unused))
121 #    else
122 #        define gmx_unused
123 #    endif
124
125 typedef float fvec[3];
126
127 // Data structures shared between OpenCL device code and OpenCL host code
128 // TODO: review, improve
129 // Replaced real by float for now, to avoid including any other header
130 typedef struct
131 {
132     float c2;
133     float c3;
134     float cpot;
135 } shift_consts_t;
136
137 /* Used with potential switching:
138  * rsw        = max(r - r_switch, 0)
139  * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
140  * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
141  * force      = force*dsw - potential*sw
142  * potential *= sw
143  */
144 typedef struct
145 {
146     float c3;
147     float c4;
148     float c5;
149 } switch_consts_t;
150
151 // Data structure shared between the OpenCL device code and OpenCL host code
152 // Must not contain OpenCL objects (buffers)
153 typedef struct cl_nbparam_params
154 {
155
156     //! type of electrostatics, takes values from #eelCu
157     int eeltype;
158     //! type of VdW impl., takes values from #evdwCu
159     int vdwtype;
160
161     //! charge multiplication factor
162     float epsfac;
163     //! Reaction-field/plain cutoff electrostatics const.
164     float c_rf;
165     //! Reaction-field electrostatics constant
166     float two_k_rf;
167     //! Ewald/PME parameter
168     float ewald_beta;
169     //! Ewald/PME correction term substracted from the direct-space potential
170     float sh_ewald;
171     //! LJ-Ewald/PME correction term added to the correction potential
172     float sh_lj_ewald;
173     //! LJ-Ewald/PME coefficient
174     float ewaldcoeff_lj;
175
176     //! Coulomb cut-off squared
177     float rcoulomb_sq;
178
179     //! VdW cut-off squared
180     float rvdw_sq;
181     //! VdW switched cut-off
182     float rvdw_switch;
183     //! Full, outer pair-list cut-off squared
184     float rlistOuter_sq;
185     //! Inner, dynamic pruned pair-list cut-off squared  XXX: this is only needed in the pruning kernels, but for now we also pass it to the nonbondeds
186     float rlistInner_sq;
187
188     //! VdW shift dispersion constants
189     shift_consts_t dispersion_shift;
190     //! VdW shift repulsion constants
191     shift_consts_t repulsion_shift;
192     //! VdW switch constants
193     switch_consts_t vdw_switch;
194
195     /* Ewald Coulomb force table data - accessed through texture memory */
196     //! table scale/spacing
197     float coulomb_tab_scale;
198 } cl_nbparam_params_t;
199
200 typedef struct
201 {
202     //! i-super-cluster
203     int sci;
204     //! Shift vector index plus possible flags
205     int shift;
206     //! Start index into cj4
207     int cj4_ind_start;
208     //! End index into cj4
209     int cj4_ind_end;
210 } nbnxn_sci_t;
211
212 typedef struct
213 {
214     //! The i-cluster interactions mask for 1 warp
215     unsigned int imask;
216     //! Index into the exclusion array for 1 warp
217     int excl_ind;
218 } nbnxn_im_ei_t;
219
220 typedef struct
221 {
222     //! The 4 j-clusters
223     int cj[4];
224     //! The i-cluster mask data       for 2 warps
225     nbnxn_im_ei_t imei[2];
226 } nbnxn_cj4_t;
227
228
229 typedef struct
230 {
231     unsigned int pair[CL_SIZE * CL_SIZE / 2]; /* Topology exclusion interaction bits for one warp,
232                                                * each unsigned has bitS for 4*8 i clusters
233                                                */
234 } nbnxn_excl_t;
235
236 /*! i-cluster interaction mask for a super-cluster with all c_nbnxnGpuNumClusterPerSupercluster bits set */
237 __constant unsigned supercl_interaction_mask = ((1U << c_nbnxnGpuNumClusterPerSupercluster) - 1U);
238
239 /*! Minimum single precision threshold for r^2 to avoid r^-12 overflow. */
240 __constant float c_nbnxnMinDistanceSquared = NBNXM_MIN_DISTANCE_SQUARED_VALUE_FLOAT;
241
242 gmx_opencl_inline void preloadCj4Generic(__local int*        sm_cjPreload,
243                                          const __global int* gm_cj,
244                                          int                 tidxi,
245                                          int                 tidxj,
246                                          bool gmx_unused iMaskCond)
247 {
248     /* Pre-load cj into shared memory */
249 #    if defined _AMD_SOURCE_ // TODO: fix by setting c_nbnxnGpuClusterpairSplit properly
250     if (tidxj == 0 & tidxi < c_nbnxnGpuJgroupSize)
251     {
252         sm_cjPreload[tidxi] = gm_cj[tidxi];
253     }
254 #    else
255     const int c_clSize                   = CL_SIZE;
256     const int c_nbnxnGpuClusterpairSplit = 2;
257     const int c_splitClSize              = c_clSize / c_nbnxnGpuClusterpairSplit;
258     if ((tidxj == 0 | tidxj == c_splitClSize) & (tidxi < c_nbnxnGpuJgroupSize))
259     {
260         sm_cjPreload[tidxi + tidxj * c_nbnxnGpuJgroupSize / c_splitClSize] = gm_cj[tidxi];
261     }
262 #    endif
263 }
264
265
266 #    if USE_SUBGROUP_PRELOAD
267 gmx_opencl_inline int preloadCj4Subgroup(const __global int* gm_cj)
268 {
269     // loads subgroup-size # of elements (8) instead of the 4 required
270     // equivalent to *cjs = *gm_cj
271     return intel_sub_group_block_read((const __global uint*)gm_cj);
272 }
273 #    endif // USE_SUBGROUP_PRELOAD
274
275 #    if USE_SUBGROUP_PRELOAD
276 typedef size_t CjType;
277 #    else
278 typedef __local int* CjType;
279 #    endif
280
281 /*! \brief Preload cj4
282  *
283  * - For AMD we load once for a wavefront of 64 threads (on 4 threads * NTHREAD_Z)
284  * - For NVIDIA once per warp (on 2x4 threads * NTHREAD_Z)
285  * - For Intel(/USE_SUBGROUP_PRELOAD) loads into private memory(/register) instead of local memory
286  *
287  * It is the caller's responsibility to make sure that data is consumed only when
288  * it's ready. This function does not call a barrier.
289  */
290 gmx_opencl_inline void preloadCj4(CjType gmx_unused* cjs,
291                                   const __global int gmx_unused* gm_cj,
292                                   int gmx_unused tidxi,
293                                   int gmx_unused tidxj,
294                                   bool gmx_unused iMaskCond)
295 {
296 #    if USE_SUBGROUP_PRELOAD
297     *cjs = preloadCj4Subgroup(gm_cj);
298 #    elif USE_CJ_PREFETCH
299     preloadCj4Generic(*cjs, gm_cj, tidxi, tidxj, iMaskCond);
300 #    else
301     // nothing to do
302 #    endif
303 }
304
305 gmx_opencl_inline int loadCjPreload(__local int* sm_cjPreload, int jm, int gmx_unused tidxi, int gmx_unused tidxj)
306 {
307 #    if defined _AMD_SOURCE_
308     int         warpLoadOffset = 0; // TODO: fix by setting c_nbnxnGpuClusterpairSplit properly
309 #    else
310     const int c_clSize = CL_SIZE;
311     const int c_nbnxnGpuClusterpairSplit = 2;
312     const int c_splitClSize = c_clSize / c_nbnxnGpuClusterpairSplit;
313     int warpLoadOffset = (tidxj & c_splitClSize) * c_nbnxnGpuJgroupSize / c_splitClSize;
314 #    endif
315     return sm_cjPreload[jm + warpLoadOffset];
316 }
317
318 /* \brief Load a cj given a jm index.
319  *
320  * If cj4 preloading is enabled, it loads from the local memory, otherwise from global.
321  */
322 gmx_opencl_inline int
323 loadCj(CjType cjs, const __global int gmx_unused* gm_cj, int jm, int gmx_unused tidxi, int gmx_unused tidxj)
324 {
325 #    if USE_SUBGROUP_PRELOAD
326     return sub_group_broadcast(cjs, jm);
327 #    elif USE_CJ_PREFETCH
328     return loadCjPreload(cjs, jm, tidxi, tidxj);
329 #    else
330     return gm_cj[jm];
331 #    endif
332 }
333
334 /*! Convert LJ sigma,epsilon parameters to C6,C12. */
335 gmx_opencl_inline float2 convert_sigma_epsilon_to_c6_c12(const float sigma, const float epsilon)
336 {
337     const float  sigma2 = sigma * sigma;
338     const float  sigma6 = sigma2 * sigma2 * sigma2;
339     const float  c6     = epsilon * sigma6;
340     const float2 c6c12  = (float2)(c6,           /* c6 */
341                                   c6 * sigma6); /* c12 */
342     return c6c12;
343 }
344
345
346 /*! Apply force switch,  force + energy version. */
347 gmx_opencl_inline void calculate_force_switch_F(const cl_nbparam_params_t* nbparam,
348                                                 float                      c6,
349                                                 float                      c12,
350                                                 float                      inv_r,
351                                                 float                      r2,
352                                                 float*                     F_invr)
353 {
354     /* force switch constants */
355     const float disp_shift_V2 = nbparam->dispersion_shift.c2;
356     const float disp_shift_V3 = nbparam->dispersion_shift.c3;
357     const float repu_shift_V2 = nbparam->repulsion_shift.c2;
358     const float repu_shift_V3 = nbparam->repulsion_shift.c3;
359
360     const float r        = r2 * inv_r;
361     float       r_switch = r - nbparam->rvdw_switch;
362     r_switch             = r_switch >= 0.0F ? r_switch : 0.0F;
363
364     *F_invr += -c6 * (disp_shift_V2 + disp_shift_V3 * r_switch) * r_switch * r_switch * inv_r
365                + c12 * (repu_shift_V2 + repu_shift_V3 * r_switch) * r_switch * r_switch * inv_r;
366 }
367
368 /*! Apply force switch, force-only version. */
369 gmx_opencl_inline void calculate_force_switch_F_E(const cl_nbparam_params_t* nbparam,
370                                                   float                      c6,
371                                                   float                      c12,
372                                                   float                      inv_r,
373                                                   float                      r2,
374                                                   float*                     F_invr,
375                                                   float*                     E_lj)
376 {
377     /* force switch constants */
378     const float disp_shift_V2 = nbparam->dispersion_shift.c2;
379     const float disp_shift_V3 = nbparam->dispersion_shift.c3;
380     const float repu_shift_V2 = nbparam->repulsion_shift.c2;
381     const float repu_shift_V3 = nbparam->repulsion_shift.c3;
382
383     const float disp_shift_F2 = nbparam->dispersion_shift.c2 / 3;
384     const float disp_shift_F3 = nbparam->dispersion_shift.c3 / 4;
385     const float repu_shift_F2 = nbparam->repulsion_shift.c2 / 3;
386     const float repu_shift_F3 = nbparam->repulsion_shift.c3 / 4;
387
388     const float r        = r2 * inv_r;
389     float       r_switch = r - nbparam->rvdw_switch;
390     r_switch             = r_switch >= 0.0F ? r_switch : 0.0F;
391
392     *F_invr += -c6 * (disp_shift_V2 + disp_shift_V3 * r_switch) * r_switch * r_switch * inv_r
393                + c12 * (repu_shift_V2 + repu_shift_V3 * r_switch) * r_switch * r_switch * inv_r;
394     *E_lj += c6 * (disp_shift_F2 + disp_shift_F3 * r_switch) * r_switch * r_switch * r_switch
395              - c12 * (repu_shift_F2 + repu_shift_F3 * r_switch) * r_switch * r_switch * r_switch;
396 }
397
398 /*! Apply potential switch, force-only version. */
399 gmx_opencl_inline void calculate_potential_switch_F(const cl_nbparam_params_t* nbparam,
400                                                     float                      inv_r,
401                                                     float                      r2,
402                                                     float*                     F_invr,
403                                                     const float*               E_lj)
404 {
405     /* potential switch constants */
406     const float switch_V3 = nbparam->vdw_switch.c3;
407     const float switch_V4 = nbparam->vdw_switch.c4;
408     const float switch_V5 = nbparam->vdw_switch.c5;
409     const float switch_F2 = nbparam->vdw_switch.c3;
410     const float switch_F3 = nbparam->vdw_switch.c4;
411     const float switch_F4 = nbparam->vdw_switch.c5;
412
413     const float r        = r2 * inv_r;
414     const float r_switch = r - nbparam->rvdw_switch;
415
416     /* Unlike in the F+E kernel, conditional is faster here */
417     if (r_switch > 0.0F)
418     {
419         const float sw = 1.0F
420                          + (switch_V3 + (switch_V4 + switch_V5 * r_switch) * r_switch) * r_switch
421                                    * r_switch * r_switch;
422         const float dsw = (switch_F2 + (switch_F3 + switch_F4 * r_switch) * r_switch) * r_switch * r_switch;
423
424         *F_invr = (*F_invr) * sw - inv_r * (*E_lj) * dsw;
425     }
426 }
427
428 /*! Apply potential switch, force + energy version. */
429 gmx_opencl_inline void calculate_potential_switch_F_E(const cl_nbparam_params_t* nbparam,
430                                                       float                      inv_r,
431                                                       float                      r2,
432                                                       float*                     F_invr,
433                                                       float*                     E_lj)
434 {
435     /* potential switch constants */
436     const float switch_V3 = nbparam->vdw_switch.c3;
437     const float switch_V4 = nbparam->vdw_switch.c4;
438     const float switch_V5 = nbparam->vdw_switch.c5;
439     const float switch_F2 = nbparam->vdw_switch.c3;
440     const float switch_F3 = nbparam->vdw_switch.c4;
441     const float switch_F4 = nbparam->vdw_switch.c5;
442
443     const float r        = r2 * inv_r;
444     float       r_switch = r - nbparam->rvdw_switch;
445     r_switch             = r_switch >= 0.0F ? r_switch : 0.0F;
446
447     /* Unlike in the F-only kernel, masking is faster here */
448     const float sw =
449             1.0F + (switch_V3 + (switch_V4 + switch_V5 * r_switch) * r_switch) * r_switch * r_switch * r_switch;
450     const float dsw = (switch_F2 + (switch_F3 + switch_F4 * r_switch) * r_switch) * r_switch * r_switch;
451
452     *F_invr = (*F_invr) * sw - inv_r * (*E_lj) * dsw;
453     *E_lj *= sw;
454 }
455
456 /*! Calculate LJ-PME grid force contribution with
457  *  geometric combination rule.
458  */
459 gmx_opencl_inline void calculate_lj_ewald_comb_geom_F(__constant const float* nbfp_comb_climg2d,
460                                                       int                     typei,
461                                                       int                     typej,
462                                                       float                   r2,
463                                                       float                   inv_r2,
464                                                       float                   lje_coeff2,
465                                                       float                   lje_coeff6_6,
466                                                       float*                  F_invr)
467 {
468     const float c6grid = nbfp_comb_climg2d[2 * typei] * nbfp_comb_climg2d[2 * typej];
469
470     /* Recalculate inv_r6 without exclusion mask */
471     const float inv_r6_nm = inv_r2 * inv_r2 * inv_r2;
472     const float cr2       = lje_coeff2 * r2;
473     const float expmcr2   = exp(-cr2);
474     const float poly      = 1.0F + cr2 + HALF_F * cr2 * cr2;
475
476     /* Subtract the grid force from the total LJ force */
477     *F_invr += c6grid * (inv_r6_nm - expmcr2 * (inv_r6_nm * poly + lje_coeff6_6)) * inv_r2;
478 }
479
480 /*! Calculate LJ-PME grid force + energy contribution with
481  *  geometric combination rule.
482  */
483 gmx_opencl_inline void calculate_lj_ewald_comb_geom_F_E(__constant const float* nbfp_comb_climg2d,
484                                                         const cl_nbparam_params_t* nbparam,
485                                                         int                        typei,
486                                                         int                        typej,
487                                                         float                      r2,
488                                                         float                      inv_r2,
489                                                         float                      lje_coeff2,
490                                                         float                      lje_coeff6_6,
491                                                         float                      int_bit,
492                                                         float*                     F_invr,
493                                                         float*                     E_lj)
494 {
495     const float c6grid = nbfp_comb_climg2d[2 * typei] * nbfp_comb_climg2d[2 * typej];
496
497     /* Recalculate inv_r6 without exclusion mask */
498     const float inv_r6_nm = inv_r2 * inv_r2 * inv_r2;
499     const float cr2       = lje_coeff2 * r2;
500     const float expmcr2   = exp(-cr2);
501     const float poly      = 1.0F + cr2 + HALF_F * cr2 * cr2;
502
503     /* Subtract the grid force from the total LJ force */
504     *F_invr += c6grid * (inv_r6_nm - expmcr2 * (inv_r6_nm * poly + lje_coeff6_6)) * inv_r2;
505
506     /* Shift should be applied only to real LJ pairs */
507     const float sh_mask = nbparam->sh_lj_ewald * int_bit;
508     *E_lj += ONE_SIXTH_F * c6grid * (inv_r6_nm * (1.0F - expmcr2 * poly) + sh_mask);
509 }
510
511 /*! Calculate LJ-PME grid force + energy contribution (if E_lj != NULL) with
512  *  Lorentz-Berthelot combination rule.
513  *  We use a single F+E kernel with conditional because the performance impact
514  *  of this is pretty small and LB on the CPU is anyway very slow.
515  */
516 gmx_opencl_inline void calculate_lj_ewald_comb_LB_F_E(__constant const float*    nbfp_comb_climg2d,
517                                                       const cl_nbparam_params_t* nbparam,
518                                                       int                        typei,
519                                                       int                        typej,
520                                                       float                      r2,
521                                                       float                      inv_r2,
522                                                       float                      lje_coeff2,
523                                                       float                      lje_coeff6_6,
524                                                       float                      int_bit,
525                                                       bool                       with_E_lj,
526                                                       float*                     F_invr,
527                                                       float*                     E_lj)
528 {
529     /* sigma and epsilon are scaled to give 6*C6 */
530     const float sigma = nbfp_comb_climg2d[2 * typei] + nbfp_comb_climg2d[2 * typej];
531
532     const float epsilon = nbfp_comb_climg2d[2 * typei + 1] * nbfp_comb_climg2d[2 * typej + 1];
533
534     const float sigma2 = sigma * sigma;
535     const float c6grid = epsilon * sigma2 * sigma2 * sigma2;
536
537     /* Recalculate inv_r6 without exclusion mask */
538     const float inv_r6_nm = inv_r2 * inv_r2 * inv_r2;
539     const float cr2       = lje_coeff2 * r2;
540     const float expmcr2   = exp(-cr2);
541     const float poly      = 1.0F + cr2 + HALF_F * cr2 * cr2;
542
543     /* Subtract the grid force from the total LJ force */
544     *F_invr += c6grid * (inv_r6_nm - expmcr2 * (inv_r6_nm * poly + lje_coeff6_6)) * inv_r2;
545
546     if (with_E_lj)
547     {
548
549         /* Shift should be applied only to real LJ pairs */
550         const float sh_mask = nbparam->sh_lj_ewald * int_bit;
551         *E_lj += ONE_SIXTH_F * c6grid * (inv_r6_nm * (1.0F - expmcr2 * poly) + sh_mask);
552     }
553 }
554
555 /*! Interpolate Ewald coulomb force using the table through the tex_nbfp texture.
556  *  Original idea: from the OpenMM project
557  */
558 gmx_opencl_inline float interpolate_coulomb_force_r(__constant const float* coulomb_tab_climg2d,
559                                                     float                   r,
560                                                     float                   scale)
561 {
562     float normalized = scale * r;
563     int   index      = (int)normalized;
564     float fract2     = normalized - (float)index;
565     float fract1     = 1.0F - fract2;
566
567     return fract1 * coulomb_tab_climg2d[index] + fract2 * coulomb_tab_climg2d[index + 1];
568 }
569
570 /*! Calculate analytical Ewald correction term. */
571 gmx_opencl_inline float pmecorrF(float z2)
572 {
573     const float FN6 = -1.7357322914161492954e-8F;
574     const float FN5 = 1.4703624142580877519e-6F;
575     const float FN4 = -0.000053401640219807709149F;
576     const float FN3 = 0.0010054721316683106153F;
577     const float FN2 = -0.019278317264888380590F;
578     const float FN1 = 0.069670166153766424023F;
579     const float FN0 = -0.75225204789749321333F;
580
581     const float FD4 = 0.0011193462567257629232F;
582     const float FD3 = 0.014866955030185295499F;
583     const float FD2 = 0.11583842382862377919F;
584     const float FD1 = 0.50736591960530292870F;
585     const float FD0 = 1.0F;
586
587     const float z4 = z2 * z2;
588
589     float polyFD0 = FD4 * z4 + FD2;
590     float polyFD1 = FD3 * z4 + FD1;
591     polyFD0       = polyFD0 * z4 + FD0;
592     polyFD0       = polyFD1 * z2 + polyFD0;
593
594     polyFD0 = 1.0F / polyFD0;
595
596     float polyFN0 = FN6 * z4 + FN4;
597     float polyFN1 = FN5 * z4 + FN3;
598     polyFN0       = polyFN0 * z4 + FN2;
599     polyFN1       = polyFN1 * z4 + FN1;
600     polyFN0       = polyFN0 * z4 + FN0;
601     polyFN0       = polyFN1 * z2 + polyFN0;
602
603     return polyFN0 * polyFD0;
604 }
605
606 #    if REDUCE_SHUFFLE
607 gmx_opencl_inline void
608 reduce_force_j_shfl(float3 fin, __global float* fout, int gmx_unused tidxi, int gmx_unused tidxj, int aidx)
609 {
610     /* Only does reduction over 4 elements in cluster. Needs to be changed
611      * for CL_SIZE>4. See CUDA code for required code */
612     fin.x += intel_sub_group_shuffle_down(fin.x, fin.x, 1);
613     fin.y += intel_sub_group_shuffle_up(fin.y, fin.y, 1);
614     fin.z += intel_sub_group_shuffle_down(fin.z, fin.z, 1);
615     if ((tidxi & 1) == 1)
616     {
617         fin.x = fin.y;
618     }
619     fin.x += intel_sub_group_shuffle_down(fin.x, fin.x, 2);
620     fin.z += intel_sub_group_shuffle_up(fin.z, fin.z, 2);
621     if (tidxi == 2)
622     {
623         fin.x = fin.z;
624     }
625     if (tidxi < 3)
626     {
627         atomicAdd_g_f(&fout[3 * aidx + tidxi], fin.x);
628     }
629 }
630 #    endif
631
632 gmx_opencl_inline void
633 reduce_force_j_generic(__local float* f_buf, float3 fcj_buf, __global float* fout, int tidxi, int tidxj, int aidx)
634 {
635     int tidx                      = tidxi + tidxj * CL_SIZE;
636     f_buf[tidx]                   = fcj_buf.x;
637     f_buf[FBUF_STRIDE + tidx]     = fcj_buf.y;
638     f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
639
640     /* Split the reduction between the first 3 column threads
641        Threads with column id 0 will do the reduction for (float3).x components
642        Threads with column id 1 will do the reduction for (float3).y components
643        Threads with column id 2 will do the reduction for (float3).z components.
644        The reduction is performed for each line tidxj of f_buf. */
645     if (tidxi < 3)
646     {
647         float f = 0.0F;
648         for (int j = tidxj * CL_SIZE; j < (tidxj + 1) * CL_SIZE; j++)
649         {
650             f += f_buf[FBUF_STRIDE * tidxi + j];
651         }
652
653         atomicAdd_g_f(&fout[3 * aidx + tidxi], f);
654     }
655 }
656
657 /*! Final j-force reduction
658  */
659 gmx_opencl_inline void reduce_force_j(__local float gmx_unused* f_buf,
660                                       float3                    fcj_buf,
661                                       __global float*           fout,
662                                       int                       tidxi,
663                                       int                       tidxj,
664                                       int                       aidx)
665 {
666 #    if REDUCE_SHUFFLE
667     reduce_force_j_shfl(fcj_buf, fout, tidxi, tidxj, aidx);
668 #    else
669     reduce_force_j_generic(f_buf, fcj_buf, fout, tidxi, tidxj, aidx);
670 #    endif
671 }
672
673 #    if REDUCE_SHUFFLE
674 gmx_opencl_inline void reduce_force_i_and_shift_shfl(__private fvec  fci_buf[],
675                                                      __global float* fout,
676                                                      bool            bCalcFshift,
677                                                      int             tidxi,
678                                                      int             tidxj,
679                                                      int             sci,
680                                                      int             shift,
681                                                      __global float* fshift)
682 {
683     /* Only does reduction over 4 elements in cluster (2 per warp). Needs to be changed
684      * for CL_SIZE>4.*/
685     float2 fshift_buf = 0;
686     for (int ci_offset = 0; ci_offset < c_nbnxnGpuNumClusterPerSupercluster; ci_offset++)
687     {
688         int    aidx = (sci * c_nbnxnGpuNumClusterPerSupercluster + ci_offset) * CL_SIZE + tidxi;
689         float3 fin  = (float3) (fci_buf[ci_offset][0], fci_buf[ci_offset][1], fci_buf[ci_offset][2]);
690         fin.x += intel_sub_group_shuffle_down(fin.x, fin.x, CL_SIZE);
691         fin.y += intel_sub_group_shuffle_up(fin.y, fin.y, CL_SIZE);
692         fin.z += intel_sub_group_shuffle_down(fin.z, fin.z, CL_SIZE);
693
694         if (tidxj & 1)
695         {
696             fin.x = fin.y;
697         }
698         /* Threads 0,1 and 2,3 increment x,y for their warp */
699         atomicAdd_g_f(&fout[3 * aidx + (tidxj & 1)], fin.x);
700         if (bCalcFshift)
701         {
702             fshift_buf[0] += fin.x;
703         }
704         /* Threads 0 and 2 increment z for their warp */
705         if ((tidxj & 1) == 0)
706         {
707             atomicAdd_g_f(&fout[3 * aidx + 2], fin.z);
708             if (bCalcFshift)
709             {
710                 fshift_buf[1] += fin.z;
711             }
712         }
713     }
714     /* add up local shift forces into global mem */
715     if (bCalcFshift)
716     {
717         // Threads 0,1 and 2,3 update x,y
718         atomicAdd_g_f(&(fshift[3 * shift + (tidxj & 1)]), fshift_buf[0]);
719         // Threads 0 and 2 update z
720         if ((tidxj & 1) == 0)
721         {
722             atomicAdd_g_f(&(fshift[3 * shift + 2]), fshift_buf[1]);
723         }
724     }
725 }
726 #    endif
727
728 /*! Final i-force reduction; this implementation works only with power of two
729  *  array sizes.
730  */
731 gmx_opencl_inline void reduce_force_i_and_shift_pow2(volatile __local float* f_buf,
732                                                      __private fvec          fci_buf[],
733                                                      __global float*         fout,
734                                                      bool                    bCalcFshift,
735                                                      int                     tidxi,
736                                                      int                     tidxj,
737                                                      int                     sci,
738                                                      int                     shift,
739                                                      __global float*         fshift)
740 {
741     float fshift_buf = 0;
742     for (int ci_offset = 0; ci_offset < c_nbnxnGpuNumClusterPerSupercluster; ci_offset++)
743     {
744         int aidx = (sci * c_nbnxnGpuNumClusterPerSupercluster + ci_offset) * CL_SIZE + tidxi;
745         int tidx = tidxi + tidxj * CL_SIZE;
746         /* store i forces in shmem */
747         f_buf[tidx]                   = fci_buf[ci_offset][0];
748         f_buf[FBUF_STRIDE + tidx]     = fci_buf[ci_offset][1];
749         f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset][2];
750         barrier(CLK_LOCAL_MEM_FENCE);
751
752         /* Reduce the initial CL_SIZE values for each i atom to half
753          * every step by using CL_SIZE * i threads.
754          * Can't just use i as loop variable because than nvcc refuses to unroll.
755          */
756         int i = CL_SIZE / 2;
757         for (int j = CL_SIZE_LOG2 - 1; j > 0; j--)
758         {
759             if (tidxj < i)
760             {
761
762                 f_buf[tidxj * CL_SIZE + tidxi] += f_buf[(tidxj + i) * CL_SIZE + tidxi];
763                 f_buf[FBUF_STRIDE + tidxj * CL_SIZE + tidxi] +=
764                         f_buf[FBUF_STRIDE + (tidxj + i) * CL_SIZE + tidxi];
765                 f_buf[2 * FBUF_STRIDE + tidxj * CL_SIZE + tidxi] +=
766                         f_buf[2 * FBUF_STRIDE + (tidxj + i) * CL_SIZE + tidxi];
767             }
768             i >>= 1;
769         }
770         /* needed because
771          * a) for CL_SIZE<8: id 2 (doing z in next block) is in 2nd warp
772          * b) for all CL_SIZE a barrier is needed before f_buf is reused by next reduce_force_i call
773          * TODO: Test on Nvidia for performance difference between having the barrier here or after the atomicAdd
774          */
775         barrier(CLK_LOCAL_MEM_FENCE);
776
777         /* i == 1, last reduction step, writing to global mem */
778         /* Split the reduction between the first 3 line threads
779            Threads with line id 0 will do the reduction for (float3).x components
780            Threads with line id 1 will do the reduction for (float3).y components
781            Threads with line id 2 will do the reduction for (float3).z components. */
782         if (tidxj < 3)
783         {
784             float f = f_buf[tidxj * FBUF_STRIDE + tidxi] + f_buf[tidxj * FBUF_STRIDE + i * CL_SIZE + tidxi];
785
786             atomicAdd_g_f(&fout[3 * aidx + tidxj], f);
787
788             if (bCalcFshift)
789             {
790                 fshift_buf += f;
791             }
792         }
793     }
794     /* add up local shift forces into global mem */
795     if (bCalcFshift)
796     {
797         /* Only threads with tidxj < 3 will update fshift.
798            The threads performing the update, must be the same as the threads
799            storing the reduction result above.
800          */
801         if (tidxj < 3)
802         {
803             atomicAdd_g_f(&(fshift[3 * shift + tidxj]), fshift_buf);
804         }
805     }
806 }
807
808 /*! Final i-force reduction
809  */
810 gmx_opencl_inline void reduce_force_i_and_shift(__local float gmx_unused* f_buf,
811                                                 __private fvec            fci_buf[],
812                                                 __global float*           f,
813                                                 bool                      bCalcFshift,
814                                                 int                       tidxi,
815                                                 int                       tidxj,
816                                                 int                       sci,
817                                                 int                       shift,
818                                                 __global float*           fshift)
819 {
820 #    if REDUCE_SHUFFLE
821     reduce_force_i_and_shift_shfl(fci_buf, f, bCalcFshift, tidxi, tidxj, sci, shift, fshift);
822 #    else
823     reduce_force_i_and_shift_pow2(f_buf, fci_buf, f, bCalcFshift, tidxi, tidxj, sci, shift, fshift);
824 #    endif
825 }
826
827
828 #    if REDUCE_SHUFFLE
829 gmx_opencl_inline void reduce_energy_shfl(float                    E_lj,
830                                           float                    E_el,
831                                           volatile __global float* e_lj,
832                                           volatile __global float* e_el,
833                                           int                      tidx)
834 {
835     E_lj = sub_group_reduce_add(E_lj);
836     E_el = sub_group_reduce_add(E_el);
837     /* Should be get_sub_group_local_id()==0. Doesn't work with Intel Classic driver.
838      * To make it possible to use REDUCE_SHUFFLE with single subgroup per i-j pair
839      * (e.g. subgroup size 16 with CL_SIZE 4), either this "if" needs to be changed or
840      * the definition of WARP_SIZE (currently CL_SIZE*CL_SIZE/2) needs to be changed
841      * (by supporting c_nbnxnGpuClusterpairSplit=1). */
842     if (tidx == 0 || tidx == WARP_SIZE)
843     {
844         atomicAdd_g_f(e_lj, E_lj);
845         atomicAdd_g_f(e_el, E_el);
846     }
847 }
848 #    endif
849
850 /*! Energy reduction; this implementation works only with power of two
851  *  array sizes.
852  */
853 gmx_opencl_inline void reduce_energy_pow2(volatile __local float*  buf,
854                                           volatile __global float* e_lj,
855                                           volatile __global float* e_el,
856                                           int                      tidx)
857 {
858     int i = WARP_SIZE / 2;
859
860     /* Can't just use i as loop variable because than nvcc refuses to unroll. */
861     for (int j = WARP_SIZE_LOG2 - 1; j > 0; j--)
862     {
863         if (tidx < i)
864         {
865             buf[tidx] += buf[tidx + i];
866             buf[FBUF_STRIDE + tidx] += buf[FBUF_STRIDE + tidx + i];
867         }
868         i >>= 1;
869     }
870
871     /* last reduction step, writing to global mem */
872     if (tidx == 0)
873     {
874         float e1 = buf[tidx] + buf[tidx + i];
875         float e2 = buf[FBUF_STRIDE + tidx] + buf[FBUF_STRIDE + tidx + i];
876
877         atomicAdd_g_f(e_lj, e1);
878         atomicAdd_g_f(e_el, e2);
879     }
880 }
881
882 gmx_opencl_inline void reduce_energy(volatile __local float gmx_unused* buf,
883                                      float                              E_lj,
884                                      float                              E_el,
885                                      volatile __global float*           e_lj,
886                                      volatile __global float*           e_el,
887                                      int                                tidx)
888 {
889 #    if REDUCE_SHUFFLE
890     reduce_energy_shfl(E_lj, E_el, e_lj, e_el, tidx);
891 #    else
892     /* flush the energies to shmem and reduce them */
893     buf[tidx] = E_lj;
894     buf[FBUF_STRIDE + tidx] = E_el;
895     reduce_energy_pow2(buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
896 #    endif
897 }
898
899 gmx_opencl_inline bool gmx_sub_group_any_localmem(volatile __local int* warp_any, int widx, bool pred)
900 {
901     if (pred)
902     {
903         warp_any[widx] = 1;
904     }
905
906     bool ret = warp_any[widx];
907
908     warp_any[widx] = 0;
909
910     return ret;
911 }
912
913 //! Returns a true if predicate is true for any work item in warp
914 gmx_opencl_inline bool gmx_sub_group_any(volatile __local int gmx_unused* warp_any, int gmx_unused widx, bool pred)
915 {
916 #    if USE_SUBGROUP_ANY
917     return sub_group_any(pred);
918 #    else
919     return gmx_sub_group_any_localmem(warp_any, widx, pred);
920 #    endif
921 }
922
923 #endif /* NBNXN_OPENCL_KERNEL_UTILS_CLH */