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)))
102 __attribute__((intel_reqd_sub_group_size(WARP_SIZE))) //2*WARP_SIZE could be enabled, see comment in reduce_energy_shfl
106 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
108 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
112 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
114 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
121 cl_nbparam_params_t nbparam_params, /* IN */
122 const __global float4 *restrict xq, /* IN */
123 __global float *restrict f, /* OUT stores float3 values */
124 __global float *restrict e_lj, /* OUT */
125 __global float *restrict e_el, /* OUT */
126 __global float *restrict fshift, /* OUT stores float3 values */
128 const __global float2 *restrict lj_comb, /* IN stores float2 values */
130 const __global int *restrict atom_types, /* IN */
132 const __global float *restrict shift_vec, /* IN stores float3 values */
133 __constant float* nbfp_climg2d, /* IN */
134 __constant float* nbfp_comb_climg2d, /* IN */
135 __constant float* coulomb_tab_climg2d, /* IN */
136 const __global nbnxn_sci_t* pl_sci, /* IN */
140 __global nbnxn_cj4_t* pl_cj4, /* OUT / IN */
141 const __global nbnxn_excl_t* excl, /* IN */
142 int bCalcFshift, /* IN */
143 __local float4 *xqib /* Pointer to dyn alloc'ed shmem */
146 /* convenience variables */
147 cl_nbparam_params_t *nbparam = &nbparam_params;
149 float rcoulomb_sq = nbparam->rcoulomb_sq;
151 float2 ljcp_i, ljcp_j;
153 #ifdef VDW_CUTOFF_CHECK
154 float rvdw_sq = nbparam_params.rvdw_sq;
158 float lje_coeff2, lje_coeff6_6;
161 float two_k_rf = nbparam->two_k_rf;
164 float coulomb_tab_scale = nbparam->coulomb_tab_scale;
167 float beta2 = nbparam->ewald_beta*nbparam->ewald_beta;
168 float beta3 = nbparam->ewald_beta*nbparam->ewald_beta*nbparam->ewald_beta;
171 float rlist_sq = nbparam->rlistOuter_sq;
176 float beta = nbparam->ewald_beta;
177 float ewald_shift = nbparam->sh_ewald;
179 float c_rf = nbparam->c_rf;
180 #endif /* EL_EWALD_ANY */
181 #endif /* CALC_ENERGIES */
183 /* thread/block/warp id-s */
184 unsigned int tidxi = get_local_id(0);
185 unsigned int tidxj = get_local_id(1);
186 unsigned int tidx = get_local_id(1) * get_local_size(0) + get_local_id(0);
187 unsigned int bidx = get_group_id(0);
188 unsigned int widx = tidx / WARP_SIZE; /* warp index */
189 int sci, ci, cj, ci_offset,
191 cij4_start, cij4_end;
195 int i, jm, j4, wexcl_idx;
198 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
199 float inv_r6, c6, c12;
201 #if defined LJ_COMB_LB
202 float sigma, epsilon;
210 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
213 unsigned int wexcl, imask, mask_ji;
215 float3 xi, xj, rv, f_ij, fcj_buf;
216 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
219 /*! i-cluster interaction mask for a super-cluster with all NCL_PER_SUPERCL=8 bits set */
220 const unsigned superClInteractionMask = ((1U << NCL_PER_SUPERCL) - 1U);
222 #define LOCAL_OFFSET (xqib + NCL_PER_SUPERCL * CL_SIZE)
225 /* shmem buffer for cj, for both warps separately */
226 cjs = (__local int *)(LOCAL_OFFSET);
228 #define LOCAL_OFFSET cjs + 2 * NBNXN_GPU_JGROUP_SIZE
233 /* shmem buffer for i atom-type pre-loading */
234 __local int *atib = (__local int *)(LOCAL_OFFSET);
236 #define LOCAL_OFFSET atib + NCL_PER_SUPERCL * CL_SIZE
238 __local float2 *ljcpib = (__local float2 *)(LOCAL_OFFSET);
240 #define LOCAL_OFFSET ljcpib + NCL_PER_SUPERCL * CL_SIZE
245 /* shmem j force buffer */
246 __local float *f_buf = (__local float *)(LOCAL_OFFSET);
248 #define LOCAL_OFFSET (f_buf + CL_SIZE * CL_SIZE * 3)
250 __local float *f_buf = 0;
252 /* Local buffer used to implement __any warp vote function from CUDA.
253 volatile is used to avoid compiler optimizations for AMD builds. */
254 volatile __local uint *warp_any = (__local uint*)(LOCAL_OFFSET);
257 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
258 sci = nb_sci.sci; /* super-cluster */
259 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
260 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
262 for (i = 0; i < NCL_PER_SUPERCL; i += CL_SIZE)
264 /* Pre-load i-atom x and q into shared memory */
265 ci = sci * NCL_PER_SUPERCL + tidxj+i;
266 ai = ci * CL_SIZE + tidxi;
268 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);
269 xqbuf.w *= nbparam->epsfac;
270 xqib[(tidxj + i) * CL_SIZE + tidxi] = xqbuf;
273 /* Pre-load the i-atom types into shared memory */
274 atib[(tidxj + i) * CL_SIZE + tidxi] = atom_types[ai];
276 ljcpib[(tidxj + i) * CL_SIZE + tidxi] = lj_comb[ai];
280 /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
281 if (tidx == 0 || tidx == WARP_SIZE)
286 barrier(CLK_LOCAL_MEM_FENCE);
288 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
290 fci_buf[ci_offset] = (float3)(0.0f);
294 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
295 lje_coeff2 = nbparam->ewaldcoeff_lj*nbparam->ewaldcoeff_lj;
296 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
297 #endif /* LJ_EWALD */
304 #if defined EXCLUSION_FORCES /* Ewald or RF */
305 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
307 /* we have the diagonal: add the charge and LJ self interaction energy term */
308 for (i = 0; i < NCL_PER_SUPERCL; i++)
310 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
311 qi = xqib[i * CL_SIZE + tidxi].w;
315 E_lj += nbfp_climg2d[atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2];
316 #endif /* LJ_EWALD */
319 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
322 E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
323 #endif /* LJ_EWALD */
325 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
326 /* Correct for epsfac^2 due to adding qi^2 */
327 E_el /= nbparam->epsfac*CL_SIZE;
328 #if defined EL_RF || defined EL_CUTOFF
331 E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
333 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
335 #endif /* EXCLUSION_FORCES */
337 #endif /* CALC_ENERGIES */
339 #ifdef EXCLUSION_FORCES
340 const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
343 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
344 for (j4 = cij4_start; j4 < cij4_end; j4++)
346 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
347 imask = pl_cj4[j4].imei[widx].imask;
348 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
350 preloadCj4(cjs, pl_cj4[j4].cj, tidxi, tidxj, imask != 0u);
356 /* Unrolling this loop improves performance without pruning but
357 * with pruning it leads to slowdown.
359 * Tested with driver 1800.5
361 * TODO: check loop unrolling with NVIDIA OpenCL
363 #if !defined PRUNE_NBL && !defined _NVIDIA_SOURCE_
367 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
369 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
371 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
373 cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
374 aj = cj * CL_SIZE + tidxj;
376 /* load j atom data */
378 xj = (float3)(xqbuf.xyz);
381 typej = atom_types[aj];
383 ljcp_j = lj_comb[aj];
386 fcj_buf = (float3)(0.0f);
388 #if !defined PRUNE_NBL
391 for (i = 0; i < NCL_PER_SUPERCL; i++)
395 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
396 ai = ci * CL_SIZE + tidxi; /* i atom index */
398 /* all threads load an atom from i cluster ci into shmem! */
399 xqbuf = xqib[i * CL_SIZE + tidxi];
400 xi = (float3)(xqbuf.xyz);
402 /* distance between i and j atoms */
407 /* vote.. should code shmem serialisation, wonder what the hit will be */
413 /* If _none_ of the atoms pairs are in cutoff range,
414 the bit corresponding to the current
415 cluster-pair in imask gets set to 0. */
424 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
426 /* cutoff & exclusion check */
427 #ifdef EXCLUSION_FORCES
428 if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
430 if ((r2 < rcoulomb_sq) * int_bit)
433 /* load the rest of the i-atom parameters */
437 typei = atib[i * CL_SIZE + tidxi];
439 ljcp_i = ljcpib[i * CL_SIZE + tidxi];
441 #else /* IATYPE_SHMEM */
443 typei = atom_types[ai];
445 ljcp_i = lj_comb[ai];
447 #endif /* IATYPE_SHMEM */
448 /* LJ 6*C6 and 12*C12 */
450 c6 = nbfp_climg2d[2 * (ntypes * typei + typej)];
451 c12 = nbfp_climg2d[2 * (ntypes * typei + typej)+1];
454 c6 = ljcp_i.x * ljcp_j.x;
455 c12 = ljcp_i.y * ljcp_j.y;
457 /* LJ 2^(1/6)*sigma and 12*epsilon */
458 sigma = ljcp_i.x + ljcp_j.x;
459 epsilon = ljcp_i.y * ljcp_j.y;
460 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
461 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
463 #endif /* LJ_COMB_GEOM */
466 // Ensure distance do not become so small that r^-12 overflows
467 r2 = max(r2, NBNXN_MIN_RSQ);
470 inv_r2 = inv_r * inv_r;
471 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
472 inv_r6 = inv_r2 * inv_r2 * inv_r2;
473 #if defined EXCLUSION_FORCES
474 /* We could mask inv_r2, but with Ewald
475 * masking both inv_r6 and F_invr is faster */
477 #endif /* EXCLUSION_FORCES */
479 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
480 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
481 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot)*ONE_TWELVETH_F -
482 c6 * (inv_r6 + nbparam->dispersion_shift.cpot)*ONE_SIXTH_F);
485 #else /* ! LJ_COMB_LB || CALC_ENERGIES */
486 float sig_r = sigma*inv_r;
487 float sig_r2 = sig_r*sig_r;
488 float sig_r6 = sig_r2*sig_r2*sig_r2;
489 #if defined EXCLUSION_FORCES
491 #endif /* EXCLUSION_FORCES */
493 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
494 #endif /* ! LJ_COMB_LB || CALC_ENERGIES */
497 #ifdef LJ_FORCE_SWITCH
499 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
501 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
502 #endif /* CALC_ENERGIES */
503 #endif /* LJ_FORCE_SWITCH */
507 #ifdef LJ_EWALD_COMB_GEOM
509 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);
511 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
512 #endif /* CALC_ENERGIES */
513 #elif defined LJ_EWALD_COMB_LB
514 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
516 int_bit, true, &F_invr, &E_lj_p
519 #endif /* CALC_ENERGIES */
521 #endif /* LJ_EWALD_COMB_GEOM */
522 #endif /* LJ_EWALD */
526 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
528 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
529 #endif /* CALC_ENERGIES */
530 #endif /* LJ_POT_SWITCH */
532 #ifdef VDW_CUTOFF_CHECK
533 /* Separate VDW cut-off check to enable twin-range cut-offs
534 * (rvdw < rcoulomb <= rlist)
536 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
537 F_invr *= vdw_in_range;
539 E_lj_p *= vdw_in_range;
541 #endif /* VDW_CUTOFF_CHECK */
550 #ifdef EXCLUSION_FORCES
551 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
553 F_invr += qi * qj_f * inv_r2 * inv_r;
557 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
559 #if defined EL_EWALD_ANA
560 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
561 #elif defined EL_EWALD_TAB
562 F_invr += qi * qj_f * (int_bit*inv_r2 -
563 interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r, coulomb_tab_scale)
565 #endif /* EL_EWALD_ANA/TAB */
569 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
572 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
575 /* 1.0f - erff is faster than erfcf */
576 E_el += qi * qj_f * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
577 #endif /* EL_EWALD_ANY */
581 /* accumulate j forces in registers */
584 /* accumulate i forces in registers */
589 /* shift the mask bit by 1 */
593 /* reduce j forces */
594 reduce_force_j(f_buf, fcj_buf, f, tidxi, tidxj, aj);
598 /* Update the imask with the new one which does not contain the
599 out of range clusters anymore. */
601 pl_cj4[j4].imei[widx].imask = imask;
606 /* skip central shifts when summing shift forces */
607 if (nb_sci.shift == CENTRAL)
611 /* reduce i forces */
612 reduce_force_i_and_shift(f_buf, fci_buf, f, bCalcFshift != 0, tidxi, tidxj, sci,
613 nb_sci.shift, fshift);
616 reduce_energy(f_buf, E_lj, E_el, e_lj, e_el, tidx);
621 #undef EXCLUSION_FORCES
626 #undef USE_CJ_PREFETCH