2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2016, 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
47 #include "nbnxn_ocl_kernel_utils.clh"
49 /////////////////////////////////////////////////////////////////////////////////////////////////
51 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
52 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
56 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
57 /* Macro to control the calculation of exclusion forces in the kernel
58 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
61 * Note: convenience macro, needs to be undef-ed at the end of the file.
63 #define EXCLUSION_FORCES
66 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
67 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
71 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
72 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
77 Kernel launch parameters:
78 - #blocks = #pair lists, blockId = pair list Id
79 - #threads = CL_SIZE^2
80 - shmem = CL_SIZE^2 * sizeof(float)
82 Each thread calculates an i force-component taking one pair of i-j atoms.
84 TODO: implement 128 threads/wavefront by porting over the NTHREAD_Z/j4 loop
85 "horizontal splitting" over threads.
89 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.
90 Thus if more strings need to be appended a new macro must be written or it must be directly appended here.
92 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
95 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
97 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
101 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
103 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
110 cl_nbparam_params_t nbparam_params, /* IN */
111 const __global float4 *restrict xq, /* IN */
112 __global float *restrict f, /* stores float3 values */ /* OUT */
113 __global float *restrict e_lj, /* OUT */
114 __global float *restrict e_el, /* OUT */
115 __global float *restrict fshift, /* stores float3 values */ /* OUT */
117 const __global float2 *restrict lj_comb, /* stores float2 values */ /* IN */
119 const __global int *restrict atom_types, /* IN */
121 const __global float *restrict shift_vec, /* stores float3 values */ /* IN */
122 __constant float* nbfp_climg2d, /* IN */
123 __constant float* nbfp_comb_climg2d, /* IN */
124 __constant float* coulomb_tab_climg2d, /* IN */
125 const __global nbnxn_sci_t* pl_sci, /* IN */
129 __global nbnxn_cj4_t* pl_cj4, /* OUT / IN */
130 const __global nbnxn_excl_t* excl, /* IN */
131 int bCalcFshift, /* IN */
132 __local float4 *xqib, /* Pointer to dyn alloc'ed shmem */
133 __global float *debug_buffer /* Debug buffer, can be used with print_to_debug_buffer_f */
136 /* convenience variables */
137 cl_nbparam_params_t *nbparam = &nbparam_params;
139 float rcoulomb_sq = nbparam->rcoulomb_sq;
141 float2 ljcp_i, ljcp_j;
143 #ifdef VDW_CUTOFF_CHECK
144 float rvdw_sq = nbparam_params.rvdw_sq;
148 float lje_coeff2, lje_coeff6_6;
151 float two_k_rf = nbparam->two_k_rf;
154 float coulomb_tab_scale = nbparam->coulomb_tab_scale;
157 float beta2 = nbparam->ewald_beta*nbparam->ewald_beta;
158 float beta3 = nbparam->ewald_beta*nbparam->ewald_beta*nbparam->ewald_beta;
161 float rlist_sq = nbparam->rlist_sq;
166 float beta = nbparam->ewald_beta;
167 float ewald_shift = nbparam->sh_ewald;
169 float c_rf = nbparam->c_rf;
170 #endif /* EL_EWALD_ANY */
171 #endif /* CALC_ENERGIES */
173 /* thread/block/warp id-s */
174 unsigned int tidxi = get_local_id(0);
175 unsigned int tidxj = get_local_id(1);
176 unsigned int tidx = get_local_id(1) * get_local_size(0) + get_local_id(0);
177 unsigned int bidx = get_group_id(0);
178 unsigned int widx = tidx / WARP_SIZE; /* warp index */
179 int sci, ci, cj, ci_offset,
181 cij4_start, cij4_end;
185 int i, jm, j4, wexcl_idx;
188 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
189 float inv_r6, c6, c12;
191 #if defined LJ_COMB_LB
192 float sigma, epsilon;
200 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
203 unsigned int wexcl, imask, mask_ji;
205 float3 xi, xj, rv, f_ij, fcj_buf/*, fshift_buf*/;
207 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
210 /*! i-cluster interaction mask for a super-cluster with all NCL_PER_SUPERCL=8 bits set */
211 const unsigned superClInteractionMask = ((1U << NCL_PER_SUPERCL) - 1U);
213 /* shmem buffer for cj, for both warps separately */
214 __local int *cjs = (__local int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
215 #define LOCAL_OFFSET cjs + 2 * NBNXN_GPU_JGROUP_SIZE
219 /* shmem buffer for i atom-type pre-loading */
220 __local int *atib = (__local int *)(LOCAL_OFFSET);
222 #define LOCAL_OFFSET atib + NCL_PER_SUPERCL * CL_SIZE
224 __local float2 *ljcpib = (__local float2 *)(LOCAL_OFFSET);
226 #define LOCAL_OFFSET ljcpib + NCL_PER_SUPERCL * CL_SIZE
230 #ifndef REDUCE_SHUFFLE
231 /* shmem j force buffer */
232 __local float *f_buf = (__local float *)(LOCAL_OFFSET);
234 #define LOCAL_OFFSET f_buf + CL_SIZE * CL_SIZE * 3
236 /* Local buffer used to implement __any warp vote function from CUDA.
237 volatile is used to avoid compiler optimizations for AMD builds. */
238 volatile __local uint *warp_any = (__local uint*)(LOCAL_OFFSET);
241 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
242 sci = nb_sci.sci; /* super-cluster */
243 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
244 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
246 /* Pre-load i-atom x and q into shared memory */
247 ci = sci * NCL_PER_SUPERCL + tidxj;
248 ai = ci * CL_SIZE + tidxi;
250 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);
251 xqbuf.w *= nbparam->epsfac;
252 xqib[tidxj * CL_SIZE + tidxi] = xqbuf;
256 /* Pre-load the i-atom types into shared memory */
257 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
259 ljcpib[tidxj * CL_SIZE + tidxi] = lj_comb[ai];
262 /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
263 if(tidx==0 || tidx==32)
266 barrier(CLK_LOCAL_MEM_FENCE);
268 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
270 fci_buf[ci_offset] = (float3)(0.0f);
274 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
275 lje_coeff2 = nbparam->ewaldcoeff_lj*nbparam->ewaldcoeff_lj;
276 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
277 #endif /* LJ_EWALD */
284 #if defined EXCLUSION_FORCES /* Ewald or RF */
285 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
287 /* we have the diagonal: add the charge and LJ self interaction energy term */
288 for (i = 0; i < NCL_PER_SUPERCL; i++)
290 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
291 qi = xqib[i * CL_SIZE + tidxi].w;
295 E_lj += nbfp_climg2d[atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2];
296 #endif /* LJ_EWALD */
299 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
302 E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
303 #endif /* LJ_EWALD */
305 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
306 /* Correct for epsfac^2 due to adding qi^2 */
307 E_el /= nbparam->epsfac*CL_SIZE;
308 #if defined EL_RF || defined EL_CUTOFF
311 E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
313 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
315 #endif /* EXCLUSION_FORCES */
317 #endif /* CALC_ENERGIES */
319 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
320 for (j4 = cij4_start; j4 < cij4_end; j4++)
322 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
323 imask = pl_cj4[j4].imei[widx].imask;
324 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
330 /* Pre-load cj into shared memory on both warps separately */
331 if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
333 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
336 /* Unrolling this loop improves performance without pruning but
337 * with pruning it leads to slowdown.
339 * Tested with driver 1800.5
341 #if !defined PRUNE_NBL
345 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
347 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
349 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
351 cj = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
352 aj = cj * CL_SIZE + tidxj;
354 /* load j atom data */
356 xj = (float3)(xqbuf.xyz);
359 typej = atom_types[aj];
361 ljcp_j = lj_comb[aj];
364 fcj_buf = (float3)(0.0f);
366 #if !defined PRUNE_NBL
369 for (i = 0; i < NCL_PER_SUPERCL; i++)
373 ci_offset = i; /* i force buffer offset */
375 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
376 ai = ci * CL_SIZE + tidxi; /* i atom index */
378 /* all threads load an atom from i cluster ci into shmem! */
379 xqbuf = xqib[i * CL_SIZE + tidxi];
380 xi = (float3)(xqbuf.xyz);
382 /* distance between i and j atoms */
387 /* vote.. should code shmem serialisation, wonder what the hit will be */
391 /* If _none_ of the atoms pairs are in cutoff range,
392 the bit corresponding to the current
393 cluster-pair in imask gets set to 0. */
400 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
402 /* cutoff & exclusion check */
403 #ifdef EXCLUSION_FORCES
404 if (r2 < rcoulomb_sq *
405 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
407 if (r2 < rcoulomb_sq * int_bit)
410 /* load the rest of the i-atom parameters */
414 typei = atib[i * CL_SIZE + tidxi];
416 ljcp_i = ljcpib[i * CL_SIZE + tidxi];
418 #else /* IATYPE_SHMEM */
420 typei = atom_types[ai];
422 ljcp_i = lj_comb[ai];
424 #endif /* IATYPE_SHMEM */
425 /* LJ 6*C6 and 12*C12 */
427 c6 = nbfp_climg2d[2 * (ntypes * typei + typej)];
428 c12 = nbfp_climg2d[2 * (ntypes * typei + typej)+1];
431 c6 = ljcp_i.x * ljcp_j.x;
432 c12 = ljcp_i.y * ljcp_j.y;
434 /* LJ 2^(1/6)*sigma and 12*epsilon */
435 sigma = ljcp_i.x + ljcp_j.x;
436 epsilon = ljcp_i.y * ljcp_j.y;
437 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
438 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
440 #endif /* LJ_COMB_GEOM */
443 // Ensure distance do not become so small that r^-12 overflows
444 r2 = max(r2,NBNXN_MIN_RSQ);
447 inv_r2 = inv_r * inv_r;
448 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
449 inv_r6 = inv_r2 * inv_r2 * inv_r2;
450 #if defined EXCLUSION_FORCES
451 /* We could mask inv_r2, but with Ewald
452 * masking both inv_r6 and F_invr is faster */
454 #endif /* EXCLUSION_FORCES */
456 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
457 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
458 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot)*ONE_TWELVETH_F -
459 c6 * (inv_r6 + nbparam->dispersion_shift.cpot)*ONE_SIXTH_F);
462 #else /* ! LJ_COMB_LB || CALC_ENERGIES */
463 float sig_r = sigma*inv_r;
464 float sig_r2 = sig_r*sig_r;
465 float sig_r6 = sig_r2*sig_r2*sig_r2;
466 #if defined EXCLUSION_FORCES
468 #endif /* EXCLUSION_FORCES */
470 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
471 #endif /* ! LJ_COMB_LB || CALC_ENERGIES */
474 #ifdef LJ_FORCE_SWITCH
476 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
478 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
479 #endif /* CALC_ENERGIES */
480 #endif /* LJ_FORCE_SWITCH */
484 #ifdef LJ_EWALD_COMB_GEOM
486 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);
488 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
489 #endif /* CALC_ENERGIES */
490 #elif defined LJ_EWALD_COMB_LB
491 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
493 int_bit, true, &F_invr, &E_lj_p
496 #endif /* CALC_ENERGIES */
498 #endif /* LJ_EWALD_COMB_GEOM */
499 #endif /* LJ_EWALD */
503 calculate_potential_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
505 calculate_potential_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
506 #endif /* CALC_ENERGIES */
507 #endif /* LJ_POT_SWITCH */
509 #ifdef VDW_CUTOFF_CHECK
510 /* Separate VDW cut-off check to enable twin-range cut-offs
511 * (rvdw < rcoulomb <= rlist)
513 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
514 F_invr *= vdw_in_range;
516 E_lj_p *= vdw_in_range;
518 #endif /* VDW_CUTOFF_CHECK */
527 #ifdef EXCLUSION_FORCES
528 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
530 F_invr += qi * qj_f * inv_r2 * inv_r;
534 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
536 #if defined EL_EWALD_ANA
537 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
538 #elif defined EL_EWALD_TAB
539 F_invr += qi * qj_f * (int_bit*inv_r2 -
541 interpolate_coulomb_force_r(nbparam->coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
543 interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r, coulomb_tab_scale)
544 #endif /* USE_TEXOBJ */
546 #endif /* EL_EWALD_ANA/TAB */
550 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
553 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
556 /* 1.0f - erff is faster than erfcf */
557 E_el += qi * qj_f * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
558 #endif /* EL_EWALD_ANY */
562 /* accumulate j forces in registers */
565 /* accumulate i forces in registers */
566 fci_buf[ci_offset] += f_ij;
570 /* shift the mask bit by 1 */
574 /* reduce j forces */
576 /* store j forces in shmem */
577 f_buf[ tidx] = fcj_buf.x;
578 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
579 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
581 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
585 /* Update the imask with the new one which does not contain the
586 out of range clusters anymore. */
588 pl_cj4[j4].imei[widx].imask = imask;
593 /* skip central shifts when summing shift forces */
594 if (nb_sci.shift == CENTRAL)
601 /* reduce i forces */
602 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
604 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
606 f_buf[ tidx] = fci_buf[ci_offset].x;
607 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
608 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
609 barrier(CLK_LOCAL_MEM_FENCE);
610 reduce_force_i(f_buf, f,
611 &fshift_buf, bCalcFshift,
613 barrier(CLK_LOCAL_MEM_FENCE);
616 /* add up local shift forces into global mem */
619 /* Only threads with tidxj < 3 will update fshift.
620 The threads performing the update must be the same with the threads
621 which stored the reduction result in reduce_force_i function
624 atomicAdd_g_f(&(fshift[3 * nb_sci.shift + tidxj]), fshift_buf);
628 /* flush the energies to shmem and reduce them */
630 f_buf[FBUF_STRIDE + tidx] = E_el;
631 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
637 #undef EXCLUSION_FORCES