2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012-2018, The GROMACS development team.
5 * Copyright (c) 2019, 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.
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.
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.
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.
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.
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.
38 * \brief OpenCL non-bonded kernel.
40 * OpenCL 1.2 support is expected.
42 * \author Anca Hamuraru <anca@streamcomputing.eu>
43 * \author Szilárd Páll <pall.szilard@gmail.com>
44 * \ingroup module_nbnxm
47 /* Currently we enable CJ prefetch for AMD/NVIDIA and disable it for the "nowarp" kernel
48 * Note that this should precede the kernel_utils include.
50 #include "nbnxm_ocl_kernel_utils.clh"
52 /////////////////////////////////////////////////////////////////////////////////////////////////
54 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
55 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
59 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD \
60 || (defined EL_CUTOFF && defined CALC_ENERGIES)
61 /* Macro to control the calculation of exclusion forces in the kernel
62 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
65 * Note: convenience macro, needs to be undef-ed at the end of the file.
67 # define EXCLUSION_FORCES
70 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
71 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
75 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
76 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
81 Kernel launch parameters:
82 - #blocks = #pair lists, blockId = pair list Id
83 - #threads = CL_SIZE^2
84 - shmem = CL_SIZE^2 * sizeof(float)
86 Each thread calculates an i force-component taking one pair of i-j atoms.
88 TODO: implement 128 threads/wavefront by porting over the NTHREAD_Z/j4 loop
89 "horizontal splitting" over threads.
93 NB_KERNEL_FUNC_NAME differs from the CUDA equivalent as it is not a variadic macro due to OpenCL
94 not having a support for them, this version only takes exactly 2 arguments. Thus if more strings
95 need to be appended a new macro must be written or it must be directly appended here.
97 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
98 #ifdef cl_intel_required_subgroup_size
99 __attribute__((intel_reqd_sub_group_size(SUBGROUP_SIZE)))
102 # ifdef CALC_ENERGIES
103 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
105 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
108 # ifdef CALC_ENERGIES
109 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
111 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
118 cl_nbparam_params_t nbparam_params, /* IN */
119 const __global float4* restrict xq, /* IN */
120 __global float* restrict f, /* OUT stores float3 values */
121 __global float* restrict gmx_unused e_lj, /* OUT */
122 __global float* restrict gmx_unused e_el, /* OUT */
123 __global float* restrict fshift, /* OUT stores float3 values */
125 const __global float2* restrict lj_comb, /* IN stores float2 values */
127 const __global int* restrict atom_types, /* IN */
129 const __global float* restrict shift_vec, /* IN stores float3 values */
130 __constant const float* gmx_unused nbfp_climg2d, /* IN */
131 __constant const float* gmx_unused nbfp_comb_climg2d, /* IN */
132 __constant const float* gmx_unused coulomb_tab_climg2d, /* IN */
133 const __global nbnxn_sci_t* pl_sci, /* IN */
137 __global nbnxn_cj4_t* pl_cj4, /* OUT / IN */
138 const __global nbnxn_excl_t* excl, /* IN */
139 int bCalcFshift, /* IN */
140 __local float4* xqib /* Pointer to dyn alloc'ed shmem */
143 /* convenience variables */
144 const cl_nbparam_params_t* const nbparam = &nbparam_params;
146 const float rcoulomb_sq = nbparam->rcoulomb_sq;
147 #ifdef VDW_CUTOFF_CHECK
148 const float rvdw_sq = nbparam_params.rvdw_sq;
151 const float two_k_rf = nbparam->two_k_rf;
154 const float coulomb_tab_scale = nbparam->coulomb_tab_scale;
157 const float beta2 = nbparam->ewald_beta * nbparam->ewald_beta;
158 const float beta3 = nbparam->ewald_beta * nbparam->ewald_beta * nbparam->ewald_beta;
161 const float rlist_sq = nbparam->rlistOuter_sq;
166 const float beta = nbparam->ewald_beta;
167 const float ewald_shift = nbparam->sh_ewald;
169 const float gmx_unused c_rf = nbparam->c_rf;
170 # endif /* EL_EWALD_ANY */
171 #endif /* CALC_ENERGIES */
173 /* thread/block/warp id-s */
174 const int tidxi = get_local_id(0);
175 const int tidxj = get_local_id(1);
176 const int tidx = (int)(get_local_id(1) * get_local_size(0) + get_local_id(0));
177 const int bidx = get_group_id(0);
178 const int widx = tidx / WARP_SIZE; /* warp index */
180 /*! i-cluster interaction mask for a super-cluster with all NCL_PER_SUPERCL=8 bits set */
181 const unsigned superClInteractionMask = ((1U << NCL_PER_SUPERCL) - 1U);
183 #define LOCAL_OFFSET (xqib + NCL_PER_SUPERCL * CL_SIZE)
186 /* shmem buffer for cj, for both warps separately */
187 cjs = (__local int*)(LOCAL_OFFSET);
189 # define LOCAL_OFFSET (cjs + 2 * c_nbnxnGpuJgroupSize)
190 #endif // USE_CJ_PREFETCH
194 /* shmem buffer for i atom-type pre-loading */
195 __local int* atib = (__local int*)(LOCAL_OFFSET); //NOLINT(google-readability-casting)
197 # define LOCAL_OFFSET (atib + NCL_PER_SUPERCL * CL_SIZE)
199 __local float2* ljcpib = (__local float2*)(LOCAL_OFFSET);
201 # define LOCAL_OFFSET (ljcpib + NCL_PER_SUPERCL * CL_SIZE)
206 /* shmem j force buffer */
207 __local float* f_buf = (__local float*)(LOCAL_OFFSET);
209 # define LOCAL_OFFSET (f_buf + CL_SIZE * CL_SIZE * 3)
211 __local float* f_buf = 0;
213 #if !USE_SUBGROUP_ANY
214 /* Local buffer used to implement __any warp vote function from CUDA.
215 volatile is used to avoid compiler optimizations for AMD builds. */
216 //NOLINTNEXTLINE(google-readability-casting)
217 volatile __local int* warp_any = (__local int*)(LOCAL_OFFSET);
219 __local int gmx_unused* warp_any = 0;
223 const nbnxn_sci_t nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
224 const int sci = nb_sci.sci; /* super-cluster */
225 const int cij4_start = nb_sci.cj4_ind_start; /* first ...*/
226 const int cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
228 for (int i = 0; i < NCL_PER_SUPERCL; i += CL_SIZE)
230 /* Pre-load i-atom x and q into shared memory */
231 const int ci = sci * NCL_PER_SUPERCL + tidxj + i;
232 const int ai = ci * CL_SIZE + tidxi;
234 float4 xqbuf = xq[ai]
235 + (float4)(shift_vec[3 * nb_sci.shift], shift_vec[3 * nb_sci.shift + 1],
236 shift_vec[3 * nb_sci.shift + 2], 0.0F);
237 xqbuf.w *= nbparam->epsfac;
238 xqib[(tidxj + i) * CL_SIZE + tidxi] = xqbuf;
241 /* Pre-load the i-atom types into shared memory */
242 atib[(tidxj + i) * CL_SIZE + tidxi] = atom_types[ai];
244 ljcpib[(tidxj + i) * CL_SIZE + tidxi] = lj_comb[ai];
248 #if !USE_SUBGROUP_ANY
249 /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
250 if (tidx == 0 || tidx == WARP_SIZE)
255 barrier(CLK_LOCAL_MEM_FENCE);
257 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
258 for (int ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
260 fci_buf[ci_offset] = (float3)(0.0F);
264 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
265 const float lje_coeff2 = nbparam->ewaldcoeff_lj * nbparam->ewaldcoeff_lj;
266 const float lje_coeff6_6 = lje_coeff2 * lje_coeff2 * lje_coeff2 * ONE_SIXTH_F;
267 #endif /* LJ_EWALD */
274 # if defined EXCLUSION_FORCES /* Ewald or RF */
275 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci * NCL_PER_SUPERCL)
277 /* we have the diagonal: add the charge and LJ self interaction energy term */
278 for (int i = 0; i < NCL_PER_SUPERCL; i++)
280 # if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
281 const float qi = xqib[i * CL_SIZE + tidxi].w;
284 # if defined LJ_EWALD
285 E_lj += nbfp_climg2d[atom_types[(sci * NCL_PER_SUPERCL + i) * CL_SIZE + tidxi] * (ntypes + 1) * 2];
286 # endif /* LJ_EWALD */
289 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
292 E_lj *= HALF_F * ONE_SIXTH_F * lje_coeff6_6;
293 # endif /* LJ_EWALD */
295 # if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
296 /* Correct for epsfac^2 due to adding qi^2 */
297 E_el /= nbparam->epsfac * CL_SIZE;
298 # if defined EL_RF || defined EL_CUTOFF
299 E_el *= -HALF_F * c_rf;
301 E_el *= -beta * M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
303 # endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
305 # endif /* EXCLUSION_FORCES */
307 #endif /* CALC_ENERGIES */
309 #ifdef EXCLUSION_FORCES
310 const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
313 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
314 for (int j4 = cij4_start; j4 < cij4_end; j4++)
316 const int wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
317 unsigned int imask = pl_cj4[j4].imei[widx].imask;
318 const unsigned int wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
320 preloadCj4(&cjs, pl_cj4[j4].cj, tidxi, tidxj, imask != 0U);
326 /* Unrolling this loop improves performance without pruning but
327 * with pruning it leads to slowdown.
329 * Tested with driver 1800.5
331 * TODO: check loop unrolling with NVIDIA OpenCL
333 #if !defined PRUNE_NBL && !defined _NVIDIA_SOURCE_
336 for (int jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
338 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
340 unsigned int mask_ji = (1U << (jm * NCL_PER_SUPERCL));
342 const int cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
343 const int aj = cj * CL_SIZE + tidxj;
345 /* load j atom data */
346 const float4 xjqbuf = xq[aj];
347 const float3 xj = (float3)(xjqbuf.xyz);
348 const float qj_f = xjqbuf.w;
350 const int typej = atom_types[aj];
352 const float2 ljcp_j = lj_comb[aj];
355 float3 fcj_buf = (float3)(0.0F);
357 #if !defined PRUNE_NBL
360 for (int i = 0; i < NCL_PER_SUPERCL; i++)
364 const int gmx_unused ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
366 /* all threads load an atom from i cluster ci into shmem! */
367 const float4 xiqbuf = xqib[i * CL_SIZE + tidxi];
368 const float3 xi = (float3)(xiqbuf.xyz);
370 /* distance between i and j atoms */
371 const float3 rv = xi - xj;
372 float r2 = norm2(rv);
375 if (!gmx_sub_group_any(warp_any, widx, r2 < rlist_sq))
381 const float int_bit = (wexcl & mask_ji) ? 1.0F : 0.0F;
383 /* cutoff & exclusion check */
384 #ifdef EXCLUSION_FORCES
385 if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
387 if ((float)(r2 < rcoulomb_sq) * int_bit != 0.0F)
390 /* load the rest of the i-atom parameters */
391 const float qi = xiqbuf.w;
394 const int typei = atib[i * CL_SIZE + tidxi];
396 const float2 ljcp_i = ljcpib[i * CL_SIZE + tidxi];
398 #else /* IATYPE_SHMEM */
399 const int ai = ci * CL_SIZE + tidxi; /* i atom index */
402 const int typei = atom_types[ai];
404 const float2 ljcp_i = lj_comb[ai];
406 #endif /* IATYPE_SHMEM */
407 /* LJ 6*C6 and 12*C12 */
409 const float c6 = nbfp_climg2d[2 * (ntypes * typei + typej)];
410 const float c12 = nbfp_climg2d[2 * (ntypes * typei + typej) + 1];
413 const float c6 = ljcp_i.x * ljcp_j.x;
414 const float c12 = ljcp_i.y * ljcp_j.y;
416 /* LJ 2^(1/6)*sigma and 12*epsilon */
417 const float sigma = ljcp_i.x + ljcp_j.x;
418 const float epsilon = ljcp_i.y * ljcp_j.y;
419 # if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
421 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
423 # endif /* LJ_COMB_GEOM */
426 // Ensure distance do not become so small that r^-12 overflows
427 r2 = max(r2, NBNXN_MIN_RSQ);
429 const float inv_r = rsqrt(r2);
430 const float inv_r2 = inv_r * inv_r;
431 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
432 float inv_r6 = inv_r2 * inv_r2 * inv_r2;
433 # if defined EXCLUSION_FORCES
434 /* We could mask inv_r2, but with Ewald
435 * masking both inv_r6 and F_invr is faster */
437 # endif /* EXCLUSION_FORCES */
439 float F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
440 # if defined CALC_ENERGIES || defined LJ_POT_SWITCH
443 * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot) * ONE_TWELVETH_F
444 - c6 * (inv_r6 + nbparam->dispersion_shift.cpot) * ONE_SIXTH_F);
447 #else /* ! LJ_COMB_LB || CALC_ENERGIES */
448 const float sig_r = sigma * inv_r;
449 const float sig_r2 = sig_r * sig_r;
450 float sig_r6 = sig_r2 * sig_r2 * sig_r2;
451 # if defined EXCLUSION_FORCES
453 # endif /* EXCLUSION_FORCES */
455 float F_invr = epsilon * sig_r6 * (sig_r6 - 1.0F) * inv_r2;
456 #endif /* ! LJ_COMB_LB || CALC_ENERGIES */
459 #ifdef LJ_FORCE_SWITCH
460 # ifdef CALC_ENERGIES
461 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
463 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
464 # endif /* CALC_ENERGIES */
465 #endif /* LJ_FORCE_SWITCH */
469 # ifdef LJ_EWALD_COMB_GEOM
470 # ifdef CALC_ENERGIES
471 calculate_lj_ewald_comb_geom_F_E(
472 nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2,
473 lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
475 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2,
476 lje_coeff2, lje_coeff6_6, &F_invr);
477 # endif /* CALC_ENERGIES */
478 # elif defined LJ_EWALD_COMB_LB
479 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej,
480 r2, inv_r2, lje_coeff2, lje_coeff6_6,
481 # ifdef CALC_ENERGIES
482 int_bit, true, &F_invr, &E_lj_p
485 # endif /* CALC_ENERGIES */
487 # endif /* LJ_EWALD_COMB_GEOM */
488 #endif /* LJ_EWALD */
491 # ifdef CALC_ENERGIES
492 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
494 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
495 # endif /* CALC_ENERGIES */
496 #endif /* LJ_POT_SWITCH */
498 #ifdef VDW_CUTOFF_CHECK
499 /* Separate VDW cut-off check to enable twin-range cut-offs
500 * (rvdw < rcoulomb <= rlist)
502 const float vdw_in_range = (r2 < rvdw_sq) ? 1.0F : 0.0F;
503 F_invr *= vdw_in_range;
504 # ifdef CALC_ENERGIES
505 E_lj_p *= vdw_in_range;
507 #endif /* VDW_CUTOFF_CHECK */
516 # ifdef EXCLUSION_FORCES
517 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
519 F_invr += qi * qj_f * inv_r2 * inv_r;
523 F_invr += qi * qj_f * (int_bit * inv_r2 * inv_r - two_k_rf);
525 #if defined EL_EWALD_ANA
527 * (int_bit * inv_r2 * inv_r + pmecorrF(beta2 * r2) * beta3);
528 #elif defined EL_EWALD_TAB
531 - interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r,
534 #endif /* EL_EWALD_ANA/TAB */
538 E_el += qi * qj_f * (int_bit * inv_r - c_rf);
541 E_el += qi * qj_f * (int_bit * inv_r + HALF_F * two_k_rf * r2 - c_rf);
544 /* 1.0F - erff is faster than erfcf */
546 * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
547 # endif /* EL_EWALD_ANY */
549 const float3 f_ij = rv * F_invr;
551 /* accumulate j forces in registers */
554 /* accumulate i forces in registers */
559 /* shift the mask bit by 1 */
563 /* reduce j forces */
564 reduce_force_j(f_buf, fcj_buf, f, tidxi, tidxj, aj);
568 /* Update the imask with the new one which does not contain the
569 out of range clusters anymore. */
571 pl_cj4[j4].imei[widx].imask = imask;
576 /* skip central shifts when summing shift forces */
577 if (nb_sci.shift == CENTRAL)
581 /* reduce i forces */
582 reduce_force_i_and_shift(f_buf, fci_buf, f, bCalcFshift != 0, tidxi, tidxj, sci, nb_sci.shift, fshift);
585 reduce_energy(f_buf, E_lj, E_el, e_lj, e_el, tidx);
590 #undef EXCLUSION_FORCES
595 #undef USE_CJ_PREFETCH