2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,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.
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 non-bonded kernel.
39 * OpenCL 1.2 support is expected.
41 * \author Anca Hamuraru <anca@streamcomputing.eu>
42 * \author Szilárd Páll <pall.szilard@gmail.com>
43 * \ingroup module_mdlib
46 /* Currently we enable CJ prefetch for AMD/NVIDIA and disable it for the "nowarp" kernel
47 * Note that this should precede the kernel_utils include.
49 #if defined _AMD_SOURCE_ || defined _NVIDIA_SOURCE_
50 #define USE_CJ_PREFETCH 1
52 #define USE_CJ_PREFETCH 0
55 #include "nbnxn_ocl_kernel_utils.clh"
57 /////////////////////////////////////////////////////////////////////////////////////////////////
59 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
60 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
64 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
65 /* Macro to control the calculation of exclusion forces in the kernel
66 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
69 * Note: convenience macro, needs to be undef-ed at the end of the file.
71 #define EXCLUSION_FORCES
74 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
75 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
79 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
80 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
85 Kernel launch parameters:
86 - #blocks = #pair lists, blockId = pair list Id
87 - #threads = CL_SIZE^2
88 - shmem = CL_SIZE^2 * sizeof(float)
90 Each thread calculates an i force-component taking one pair of i-j atoms.
92 TODO: implement 128 threads/wavefront by porting over the NTHREAD_Z/j4 loop
93 "horizontal splitting" over threads.
97 NB_KERNEL_FUNC_NAME differs from the CUDA equivalent as it is not a variadic macro due to OpenCL not having a support for them, this version only takes exactly 2 arguments.
98 Thus if more strings need to be appended a new macro must be written or it must be directly appended here.
100 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
103 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
105 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
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 e_lj, /* OUT */
122 __global float *restrict 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 float* nbfp_climg2d, /* IN */
131 __constant float* nbfp_comb_climg2d, /* IN */
132 __constant float* 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 cl_nbparam_params_t *nbparam = &nbparam_params;
146 float rcoulomb_sq = nbparam->rcoulomb_sq;
148 float2 ljcp_i, ljcp_j;
150 #ifdef VDW_CUTOFF_CHECK
151 float rvdw_sq = nbparam_params.rvdw_sq;
155 float lje_coeff2, lje_coeff6_6;
158 float two_k_rf = nbparam->two_k_rf;
161 float coulomb_tab_scale = nbparam->coulomb_tab_scale;
164 float beta2 = nbparam->ewald_beta*nbparam->ewald_beta;
165 float beta3 = nbparam->ewald_beta*nbparam->ewald_beta*nbparam->ewald_beta;
168 float rlist_sq = nbparam->rlistOuter_sq;
173 float beta = nbparam->ewald_beta;
174 float ewald_shift = nbparam->sh_ewald;
176 float c_rf = nbparam->c_rf;
177 #endif /* EL_EWALD_ANY */
178 #endif /* CALC_ENERGIES */
180 /* thread/block/warp id-s */
181 unsigned int tidxi = get_local_id(0);
182 unsigned int tidxj = get_local_id(1);
183 unsigned int tidx = get_local_id(1) * get_local_size(0) + get_local_id(0);
184 unsigned int bidx = get_group_id(0);
185 unsigned int widx = tidx / WARP_SIZE; /* warp index */
186 int sci, ci, cj, ci_offset,
188 cij4_start, cij4_end;
192 int i, jm, j4, wexcl_idx;
195 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
196 float inv_r6, c6, c12;
198 #if defined LJ_COMB_LB
199 float sigma, epsilon;
207 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
210 unsigned int wexcl, imask, mask_ji;
212 float3 xi, xj, rv, f_ij, fcj_buf /*, fshift_buf*/;
214 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
217 /*! i-cluster interaction mask for a super-cluster with all NCL_PER_SUPERCL=8 bits set */
218 const unsigned superClInteractionMask = ((1U << NCL_PER_SUPERCL) - 1U);
220 #define LOCAL_OFFSET xqib + NCL_PER_SUPERCL * CL_SIZE
223 /* shmem buffer for cj, for both warps separately */
224 cjs = (__local int *)(LOCAL_OFFSET);
226 #define LOCAL_OFFSET cjs + 2 * NBNXN_GPU_JGROUP_SIZE
231 /* shmem buffer for i atom-type pre-loading */
232 __local int *atib = (__local int *)(LOCAL_OFFSET);
234 #define LOCAL_OFFSET atib + NCL_PER_SUPERCL * CL_SIZE
236 __local float2 *ljcpib = (__local float2 *)(LOCAL_OFFSET);
238 #define LOCAL_OFFSET ljcpib + NCL_PER_SUPERCL * CL_SIZE
242 #ifndef REDUCE_SHUFFLE
243 /* shmem j force buffer */
244 __local float *f_buf = (__local float *)(LOCAL_OFFSET);
246 #define LOCAL_OFFSET f_buf + CL_SIZE * CL_SIZE * 3
248 /* Local buffer used to implement __any warp vote function from CUDA.
249 volatile is used to avoid compiler optimizations for AMD builds. */
250 volatile __local uint *warp_any = (__local uint*)(LOCAL_OFFSET);
253 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
254 sci = nb_sci.sci; /* super-cluster */
255 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
256 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
258 for (i = 0; i < NCL_PER_SUPERCL; i += CL_SIZE)
260 /* Pre-load i-atom x and q into shared memory */
261 ci = sci * NCL_PER_SUPERCL + tidxj+i;
262 ai = ci * CL_SIZE + tidxi;
264 xqbuf = xq[ai] + (float4)(shift_vec[3 * nb_sci.shift], shift_vec[3 * nb_sci.shift + 1], shift_vec[3 * nb_sci.shift + 2], 0.0f);
265 xqbuf.w *= nbparam->epsfac;
266 xqib[(tidxj + i) * CL_SIZE + tidxi] = xqbuf;
269 /* Pre-load the i-atom types into shared memory */
270 atib[(tidxj + i) * CL_SIZE + tidxi] = atom_types[ai];
272 ljcpib[(tidxj + i) * CL_SIZE + tidxi] = lj_comb[ai];
276 /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
277 if (tidx == 0 || tidx == WARP_SIZE)
282 barrier(CLK_LOCAL_MEM_FENCE);
284 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
286 fci_buf[ci_offset] = (float3)(0.0f);
290 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
291 lje_coeff2 = nbparam->ewaldcoeff_lj*nbparam->ewaldcoeff_lj;
292 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
293 #endif /* LJ_EWALD */
300 #if defined EXCLUSION_FORCES /* Ewald or RF */
301 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
303 /* we have the diagonal: add the charge and LJ self interaction energy term */
304 for (i = 0; i < NCL_PER_SUPERCL; i++)
306 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
307 qi = xqib[i * CL_SIZE + tidxi].w;
311 E_lj += nbfp_climg2d[atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2];
312 #endif /* LJ_EWALD */
315 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
318 E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
319 #endif /* LJ_EWALD */
321 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
322 /* Correct for epsfac^2 due to adding qi^2 */
323 E_el /= nbparam->epsfac*CL_SIZE;
324 #if defined EL_RF || defined EL_CUTOFF
327 E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
329 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
331 #endif /* EXCLUSION_FORCES */
333 #endif /* CALC_ENERGIES */
335 #ifdef EXCLUSION_FORCES
336 const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
339 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
340 for (j4 = cij4_start; j4 < cij4_end; j4++)
342 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
343 imask = pl_cj4[j4].imei[widx].imask;
344 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
346 preloadCj4(cjs, pl_cj4[j4].cj, tidxi, tidxj, imask);
352 /* Unrolling this loop improves performance without pruning but
353 * with pruning it leads to slowdown.
355 * Tested with driver 1800.5
357 * TODO: check loop unrolling with NVIDIA OpenCL
359 #if !defined PRUNE_NBL && !defined _NVIDIA_SOURCE_
363 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
365 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
367 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
369 cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
370 aj = cj * CL_SIZE + tidxj;
372 /* load j atom data */
374 xj = (float3)(xqbuf.xyz);
377 typej = atom_types[aj];
379 ljcp_j = lj_comb[aj];
382 fcj_buf = (float3)(0.0f);
384 #if !defined PRUNE_NBL
387 for (i = 0; i < NCL_PER_SUPERCL; i++)
391 ci_offset = i; /* i force buffer offset */
393 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
394 ai = ci * CL_SIZE + tidxi; /* i atom index */
396 /* all threads load an atom from i cluster ci into shmem! */
397 xqbuf = xqib[i * CL_SIZE + tidxi];
398 xi = (float3)(xqbuf.xyz);
400 /* distance between i and j atoms */
405 /* vote.. should code shmem serialisation, wonder what the hit will be */
411 /* If _none_ of the atoms pairs are in cutoff range,
412 the bit corresponding to the current
413 cluster-pair in imask gets set to 0. */
422 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
424 /* cutoff & exclusion check */
425 #ifdef EXCLUSION_FORCES
426 if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
428 if ((r2 < rcoulomb_sq) * int_bit)
431 /* load the rest of the i-atom parameters */
435 typei = atib[i * CL_SIZE + tidxi];
437 ljcp_i = ljcpib[i * CL_SIZE + tidxi];
439 #else /* IATYPE_SHMEM */
441 typei = atom_types[ai];
443 ljcp_i = lj_comb[ai];
445 #endif /* IATYPE_SHMEM */
446 /* LJ 6*C6 and 12*C12 */
448 c6 = nbfp_climg2d[2 * (ntypes * typei + typej)];
449 c12 = nbfp_climg2d[2 * (ntypes * typei + typej)+1];
452 c6 = ljcp_i.x * ljcp_j.x;
453 c12 = ljcp_i.y * ljcp_j.y;
455 /* LJ 2^(1/6)*sigma and 12*epsilon */
456 sigma = ljcp_i.x + ljcp_j.x;
457 epsilon = ljcp_i.y * ljcp_j.y;
458 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
459 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
461 #endif /* LJ_COMB_GEOM */
464 // Ensure distance do not become so small that r^-12 overflows
465 r2 = max(r2, NBNXN_MIN_RSQ);
468 inv_r2 = inv_r * inv_r;
469 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
470 inv_r6 = inv_r2 * inv_r2 * inv_r2;
471 #if defined EXCLUSION_FORCES
472 /* We could mask inv_r2, but with Ewald
473 * masking both inv_r6 and F_invr is faster */
475 #endif /* EXCLUSION_FORCES */
477 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
478 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
479 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot)*ONE_TWELVETH_F -
480 c6 * (inv_r6 + nbparam->dispersion_shift.cpot)*ONE_SIXTH_F);
483 #else /* ! LJ_COMB_LB || CALC_ENERGIES */
484 float sig_r = sigma*inv_r;
485 float sig_r2 = sig_r*sig_r;
486 float sig_r6 = sig_r2*sig_r2*sig_r2;
487 #if defined EXCLUSION_FORCES
489 #endif /* EXCLUSION_FORCES */
491 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
492 #endif /* ! LJ_COMB_LB || CALC_ENERGIES */
495 #ifdef LJ_FORCE_SWITCH
497 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
499 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
500 #endif /* CALC_ENERGIES */
501 #endif /* LJ_FORCE_SWITCH */
505 #ifdef LJ_EWALD_COMB_GEOM
507 calculate_lj_ewald_comb_geom_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
509 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
510 #endif /* CALC_ENERGIES */
511 #elif defined LJ_EWALD_COMB_LB
512 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
514 int_bit, true, &F_invr, &E_lj_p
517 #endif /* CALC_ENERGIES */
519 #endif /* LJ_EWALD_COMB_GEOM */
520 #endif /* LJ_EWALD */
524 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
526 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
527 #endif /* CALC_ENERGIES */
528 #endif /* LJ_POT_SWITCH */
530 #ifdef VDW_CUTOFF_CHECK
531 /* Separate VDW cut-off check to enable twin-range cut-offs
532 * (rvdw < rcoulomb <= rlist)
534 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
535 F_invr *= vdw_in_range;
537 E_lj_p *= vdw_in_range;
539 #endif /* VDW_CUTOFF_CHECK */
548 #ifdef EXCLUSION_FORCES
549 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
551 F_invr += qi * qj_f * inv_r2 * inv_r;
555 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
557 #if defined EL_EWALD_ANA
558 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
559 #elif defined EL_EWALD_TAB
560 F_invr += qi * qj_f * (int_bit*inv_r2 -
561 interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r, coulomb_tab_scale)
563 #endif /* EL_EWALD_ANA/TAB */
567 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
570 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
573 /* 1.0f - erff is faster than erfcf */
574 E_el += qi * qj_f * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
575 #endif /* EL_EWALD_ANY */
579 /* accumulate j forces in registers */
582 /* accumulate i forces in registers */
583 fci_buf[ci_offset] += f_ij;
587 /* shift the mask bit by 1 */
591 /* reduce j forces */
593 /* store j forces in shmem */
594 f_buf[ tidx] = fcj_buf.x;
595 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
596 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
598 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
602 /* Update the imask with the new one which does not contain the
603 out of range clusters anymore. */
605 pl_cj4[j4].imei[widx].imask = imask;
610 /* skip central shifts when summing shift forces */
611 if (nb_sci.shift == CENTRAL)
618 /* reduce i forces */
619 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
621 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
623 f_buf[ tidx] = fci_buf[ci_offset].x;
624 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
625 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
626 barrier(CLK_LOCAL_MEM_FENCE);
627 reduce_force_i(f_buf, f,
628 &fshift_buf, bCalcFshift,
632 /* add up local shift forces into global mem */
635 /* Only threads with tidxj < 3 will update fshift.
636 The threads performing the update must be the same with the threads
637 which stored the reduction result in reduce_force_i function
641 atomicAdd_g_f(&(fshift[3 * nb_sci.shift + tidxj]), fshift_buf);
646 /* flush the energies to shmem and reduce them */
648 f_buf[FBUF_STRIDE + tidx] = E_el;
649 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
654 #undef EXCLUSION_FORCES
659 #undef USE_CJ_PREFETCH