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 for AMD GPUs.
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 #define USE_CJ_PREFETCH 1
51 #include "nbnxn_ocl_kernel_utils.clh"
53 /////////////////////////////////////////////////////////////////////////////////////////////////
55 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
56 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
60 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (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 not having a support for them, this version only takes exactly 2 arguments.
94 Thus if more strings need to be appended a new macro must be written or it must be directly appended here.
96 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
99 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
101 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
105 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
107 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
114 cl_nbparam_params_t nbparam_params, /* IN */
115 const __global float4 *restrict xq, /* IN */
116 __global float *restrict f, /* OUT stores float3 values */
117 __global float *restrict e_lj, /* OUT */
118 __global float *restrict e_el, /* OUT */
119 __global float *restrict fshift, /* OUT stores float3 values */
121 const __global float2 *restrict lj_comb, /* IN stores float2 values */
123 const __global int *restrict atom_types, /* IN */
125 const __global float *restrict shift_vec, /* IN stores float3 values */
126 __constant float* nbfp_climg2d, /* IN */
127 __constant float* nbfp_comb_climg2d, /* IN */
128 __constant float* coulomb_tab_climg2d, /* IN */
129 const __global nbnxn_sci_t* pl_sci, /* IN */
133 __global nbnxn_cj4_t* pl_cj4, /* OUT / IN */
134 const __global nbnxn_excl_t* excl, /* IN */
135 int bCalcFshift, /* IN */
136 __local float4 *xqib /* Pointer to dyn alloc'ed shmem */
139 /* convenience variables */
140 cl_nbparam_params_t *nbparam = &nbparam_params;
142 float rcoulomb_sq = nbparam->rcoulomb_sq;
144 float2 ljcp_i, ljcp_j;
146 #ifdef VDW_CUTOFF_CHECK
147 float rvdw_sq = nbparam_params.rvdw_sq;
151 float lje_coeff2, lje_coeff6_6;
154 float two_k_rf = nbparam->two_k_rf;
157 float coulomb_tab_scale = nbparam->coulomb_tab_scale;
160 float beta2 = nbparam->ewald_beta*nbparam->ewald_beta;
161 float beta3 = nbparam->ewald_beta*nbparam->ewald_beta*nbparam->ewald_beta;
164 float rlist_sq = nbparam->rlistOuter_sq;
169 float beta = nbparam->ewald_beta;
170 float ewald_shift = nbparam->sh_ewald;
172 float c_rf = nbparam->c_rf;
173 #endif /* EL_EWALD_ANY */
174 #endif /* CALC_ENERGIES */
176 /* thread/block/warp id-s */
177 unsigned int tidxi = get_local_id(0);
178 unsigned int tidxj = get_local_id(1);
179 unsigned int tidx = get_local_id(1) * get_local_size(0) + get_local_id(0);
180 unsigned int bidx = get_group_id(0);
181 unsigned int widx = tidx / WARP_SIZE; /* warp index */
182 int sci, ci, cj, ci_offset,
184 cij4_start, cij4_end;
188 int i, jm, j4, wexcl_idx;
191 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
192 float inv_r6, c6, c12;
194 #if defined LJ_COMB_LB
195 float sigma, epsilon;
203 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
206 unsigned int wexcl, imask, mask_ji;
208 float3 xi, xj, rv, f_ij, fcj_buf /*, fshift_buf*/;
210 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
213 /*! i-cluster interaction mask for a super-cluster with all NCL_PER_SUPERCL=8 bits set */
214 const unsigned superClInteractionMask = ((1U << NCL_PER_SUPERCL) - 1U);
216 #define LOCAL_OFFSET xqib + NCL_PER_SUPERCL * CL_SIZE
219 /* shmem buffer for cj, for both warps separately */
220 cjs = (__local int *)(LOCAL_OFFSET);
222 #define LOCAL_OFFSET cjs + 2 * NBNXN_GPU_JGROUP_SIZE
227 /* shmem buffer for i atom-type pre-loading */
228 __local int *atib = (__local int *)(LOCAL_OFFSET);
230 #define LOCAL_OFFSET atib + NCL_PER_SUPERCL * CL_SIZE
232 __local float2 *ljcpib = (__local float2 *)(LOCAL_OFFSET);
234 #define LOCAL_OFFSET ljcpib + NCL_PER_SUPERCL * CL_SIZE
238 #ifndef REDUCE_SHUFFLE
239 /* shmem j force buffer */
240 __local float *f_buf = (__local float *)(LOCAL_OFFSET);
242 #define LOCAL_OFFSET f_buf + CL_SIZE * CL_SIZE * 3
244 /* Local buffer used to implement __any warp vote function from CUDA.
245 volatile is used to avoid compiler optimizations for AMD builds. */
246 volatile __local uint *warp_any = (__local uint*)(LOCAL_OFFSET);
249 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
250 sci = nb_sci.sci; /* super-cluster */
251 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
252 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
254 /* Pre-load i-atom x and q into shared memory */
255 ci = sci * NCL_PER_SUPERCL + tidxj;
256 ai = ci * CL_SIZE + tidxi;
258 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);
259 xqbuf.w *= nbparam->epsfac;
260 xqib[tidxj * CL_SIZE + tidxi] = xqbuf;
264 /* Pre-load the i-atom types into shared memory */
265 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
267 ljcpib[tidxj * CL_SIZE + tidxi] = lj_comb[ai];
270 /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
271 if (tidx == 0 || tidx == 32)
276 barrier(CLK_LOCAL_MEM_FENCE);
278 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
280 fci_buf[ci_offset] = (float3)(0.0f);
284 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
285 lje_coeff2 = nbparam->ewaldcoeff_lj*nbparam->ewaldcoeff_lj;
286 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
287 #endif /* LJ_EWALD */
294 #if defined EXCLUSION_FORCES /* Ewald or RF */
295 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
297 /* we have the diagonal: add the charge and LJ self interaction energy term */
298 for (i = 0; i < NCL_PER_SUPERCL; i++)
300 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
301 qi = xqib[i * CL_SIZE + tidxi].w;
305 E_lj += nbfp_climg2d[atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2];
306 #endif /* LJ_EWALD */
309 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
312 E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
313 #endif /* LJ_EWALD */
315 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
316 /* Correct for epsfac^2 due to adding qi^2 */
317 E_el /= nbparam->epsfac*CL_SIZE;
318 #if defined EL_RF || defined EL_CUTOFF
321 E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
323 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
325 #endif /* EXCLUSION_FORCES */
327 #endif /* CALC_ENERGIES */
329 #ifdef EXCLUSION_FORCES
330 const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
333 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
334 for (j4 = cij4_start; j4 < cij4_end; j4++)
336 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
337 imask = pl_cj4[j4].imei[widx].imask;
338 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
340 preloadCj4(cjs, pl_cj4[j4].cj, tidxi, tidxj, imask);
346 /* Unrolling this loop improves performance without pruning but
347 * with pruning it leads to slowdown.
349 * Tested with driver 1800.5
351 #if !defined PRUNE_NBL
355 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
357 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
359 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
361 cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
362 aj = cj * CL_SIZE + tidxj;
364 /* load j atom data */
366 xj = (float3)(xqbuf.xyz);
369 typej = atom_types[aj];
371 ljcp_j = lj_comb[aj];
374 fcj_buf = (float3)(0.0f);
376 #if !defined PRUNE_NBL
379 for (i = 0; i < NCL_PER_SUPERCL; i++)
383 ci_offset = i; /* i force buffer offset */
385 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
386 ai = ci * CL_SIZE + tidxi; /* i atom index */
388 /* all threads load an atom from i cluster ci into shmem! */
389 xqbuf = xqib[i * CL_SIZE + tidxi];
390 xi = (float3)(xqbuf.xyz);
392 /* distance between i and j atoms */
397 /* vote.. should code shmem serialisation, wonder what the hit will be */
403 /* If _none_ of the atoms pairs are in cutoff range,
404 the bit corresponding to the current
405 cluster-pair in imask gets set to 0. */
414 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
416 /* cutoff & exclusion check */
417 #ifdef EXCLUSION_FORCES
418 if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
420 if ((r2 < rcoulomb_sq) * int_bit)
423 /* load the rest of the i-atom parameters */
427 typei = atib[i * CL_SIZE + tidxi];
429 ljcp_i = ljcpib[i * CL_SIZE + tidxi];
431 #else /* IATYPE_SHMEM */
433 typei = atom_types[ai];
435 ljcp_i = lj_comb[ai];
437 #endif /* IATYPE_SHMEM */
438 /* LJ 6*C6 and 12*C12 */
440 c6 = nbfp_climg2d[2 * (ntypes * typei + typej)];
441 c12 = nbfp_climg2d[2 * (ntypes * typei + typej)+1];
444 c6 = ljcp_i.x * ljcp_j.x;
445 c12 = ljcp_i.y * ljcp_j.y;
447 /* LJ 2^(1/6)*sigma and 12*epsilon */
448 sigma = ljcp_i.x + ljcp_j.x;
449 epsilon = ljcp_i.y * ljcp_j.y;
450 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
451 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
453 #endif /* LJ_COMB_GEOM */
456 // Ensure distance do not become so small that r^-12 overflows
457 r2 = max(r2, NBNXN_MIN_RSQ);
460 inv_r2 = inv_r * inv_r;
461 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
462 inv_r6 = inv_r2 * inv_r2 * inv_r2;
463 #if defined EXCLUSION_FORCES
464 /* We could mask inv_r2, but with Ewald
465 * masking both inv_r6 and F_invr is faster */
467 #endif /* EXCLUSION_FORCES */
469 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
470 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
471 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot)*ONE_TWELVETH_F -
472 c6 * (inv_r6 + nbparam->dispersion_shift.cpot)*ONE_SIXTH_F);
475 #else /* ! LJ_COMB_LB || CALC_ENERGIES */
476 float sig_r = sigma*inv_r;
477 float sig_r2 = sig_r*sig_r;
478 float sig_r6 = sig_r2*sig_r2*sig_r2;
479 #if defined EXCLUSION_FORCES
481 #endif /* EXCLUSION_FORCES */
483 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
484 #endif /* ! LJ_COMB_LB || CALC_ENERGIES */
487 #ifdef LJ_FORCE_SWITCH
489 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
491 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
492 #endif /* CALC_ENERGIES */
493 #endif /* LJ_FORCE_SWITCH */
497 #ifdef LJ_EWALD_COMB_GEOM
499 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);
501 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
502 #endif /* CALC_ENERGIES */
503 #elif defined LJ_EWALD_COMB_LB
504 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
506 int_bit, true, &F_invr, &E_lj_p
509 #endif /* CALC_ENERGIES */
511 #endif /* LJ_EWALD_COMB_GEOM */
512 #endif /* LJ_EWALD */
516 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
518 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
519 #endif /* CALC_ENERGIES */
520 #endif /* LJ_POT_SWITCH */
522 #ifdef VDW_CUTOFF_CHECK
523 /* Separate VDW cut-off check to enable twin-range cut-offs
524 * (rvdw < rcoulomb <= rlist)
526 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
527 F_invr *= vdw_in_range;
529 E_lj_p *= vdw_in_range;
531 #endif /* VDW_CUTOFF_CHECK */
540 #ifdef EXCLUSION_FORCES
541 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
543 F_invr += qi * qj_f * inv_r2 * inv_r;
547 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
549 #if defined EL_EWALD_ANA
550 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
551 #elif defined EL_EWALD_TAB
552 F_invr += qi * qj_f * (int_bit*inv_r2 -
553 interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r, coulomb_tab_scale)
555 #endif /* EL_EWALD_ANA/TAB */
559 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
562 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
565 /* 1.0f - erff is faster than erfcf */
566 E_el += qi * qj_f * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
567 #endif /* EL_EWALD_ANY */
571 /* accumulate j forces in registers */
574 /* accumulate i forces in registers */
575 fci_buf[ci_offset] += f_ij;
579 /* shift the mask bit by 1 */
583 /* reduce j forces */
585 /* store j forces in shmem */
586 f_buf[ tidx] = fcj_buf.x;
587 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
588 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
590 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
594 /* Update the imask with the new one which does not contain the
595 out of range clusters anymore. */
597 pl_cj4[j4].imei[widx].imask = imask;
602 /* skip central shifts when summing shift forces */
603 if (nb_sci.shift == CENTRAL)
610 /* reduce i forces */
611 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
613 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
615 f_buf[ tidx] = fci_buf[ci_offset].x;
616 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
617 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
618 barrier(CLK_LOCAL_MEM_FENCE);
619 reduce_force_i(f_buf, f,
620 &fshift_buf, bCalcFshift,
622 barrier(CLK_LOCAL_MEM_FENCE);
625 /* add up local shift forces into global mem */
628 /* Only threads with tidxj < 3 will update fshift.
629 The threads performing the update must be the same with the threads
630 which stored the reduction result in reduce_force_i function
634 atomicAdd_g_f(&(fshift[3 * nb_sci.shift + tidxj]), fshift_buf);
639 /* flush the energies to shmem and reduce them */
641 f_buf[FBUF_STRIDE + tidx] = E_el;
642 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
648 #undef EXCLUSION_FORCES
653 #undef USE_CJ_PREFETCH