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.
36 #include "nbnxn_ocl_kernel_utils.clh"
38 /////////////////////////////////////////////////////////////////////////////////////////////////
40 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
41 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
45 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
46 /* Macro to control the calculation of exclusion forces in the kernel
47 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
50 * Note: convenience macro, needs to be undef-ed at the end of the file.
52 #define EXCLUSION_FORCES
55 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
56 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
61 Kernel launch parameters:
62 - #blocks = #pair lists, blockId = pair list Id
63 - #threads = CL_SIZE^2
64 - shmem = CL_SIZE^2 * sizeof(float)
66 Each thread calculates an i force-component taking one pair of i-j atoms.
68 TODO: implement 128 threads/wavefront by porting over the NTHREAD_Z/j4 loop
69 "horizontal splitting" over threads.
73 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.
74 Thus if more strings need to be appended a new macro must be written or it must be directly appended here.
76 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
79 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
81 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
85 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
87 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
91 cl_nbparam_params_t nbparam_params, /* IN */
92 const __global float4 *restrict xq, /* IN */
93 __global float *restrict f, /* stores float3 values */ /* OUT */
94 __global float *restrict e_lj, /* OUT */
95 __global float *restrict e_el, /* OUT */
96 __global float *restrict fshift, /* stores float3 values */ /* OUT */
97 const __global int *restrict atom_types, /* IN */
98 const __global float *restrict shift_vec, /* stores float3 values */ /* IN */
99 __constant float* nbfp_climg2d, /* IN */
100 __constant float* nbfp_comb_climg2d, /* IN */
101 __constant float* coulomb_tab_climg2d, /* IN */
102 const __global nbnxn_sci_t* pl_sci, /* IN */
106 __global nbnxn_cj4_t* pl_cj4, /* OUT / IN */
107 const __global nbnxn_excl_t* excl, /* IN */
108 int bCalcFshift, /* IN */
109 __local float4 *xqib, /* Pointer to dyn alloc'ed shmem */
110 __global float *debug_buffer /* Debug buffer, can be used with print_to_debug_buffer_f */
113 /* convenience variables */
114 cl_nbparam_params_t *nbparam = &nbparam_params;
116 float rcoulomb_sq = nbparam->rcoulomb_sq;
118 #ifdef VDW_CUTOFF_CHECK
119 float rvdw_sq = nbparam_params.rvdw_sq;
123 float lje_coeff2, lje_coeff6_6;
126 float two_k_rf = nbparam->two_k_rf;
129 float coulomb_tab_scale = nbparam->coulomb_tab_scale;
132 float beta2 = nbparam->ewald_beta*nbparam->ewald_beta;
133 float beta3 = nbparam->ewald_beta*nbparam->ewald_beta*nbparam->ewald_beta;
136 float rlist_sq = nbparam->rlist_sq;
141 float beta = nbparam->ewald_beta;
142 float ewald_shift = nbparam->sh_ewald;
144 float c_rf = nbparam->c_rf;
145 #endif /* EL_EWALD_ANY */
146 #endif /* CALC_ENERGIES */
148 /* thread/block/warp id-s */
149 unsigned int tidxi = get_local_id(0);
150 unsigned int tidxj = get_local_id(1);
151 unsigned int tidx = get_local_id(1) * get_local_size(0) + get_local_id(0);
152 unsigned int bidx = get_group_id(0);
153 unsigned int widx = tidx / WARP_SIZE; /* warp index */
154 int sci, ci, cj, ci_offset,
156 cij4_start, cij4_end,
158 i, jm, j4, wexcl_idx;
160 r2, inv_r, inv_r2, inv_r6,
168 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
171 unsigned int wexcl, imask, mask_ji;
173 float3 xi, xj, rv, f_ij, fcj_buf/*, fshift_buf*/;
175 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
178 /*! i-cluster interaction mask for a super-cluster with all NCL_PER_SUPERCL=8 bits set */
179 const unsigned superClInteractionMask = ((1U << NCL_PER_SUPERCL) - 1U);
181 /* shmem buffer for cj, for both warps separately */
182 __local int *cjs = (__local int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
183 #define LOCAL_OFFSET cjs + 2 * NBNXN_GPU_JGROUP_SIZE
185 #ifdef IATYPE_SHMEM //Should not be defined!
186 /* shmem buffer for i atom-type pre-loading */
187 __local int *atib = (__local int *)(LOCAL_OFFSET);
189 #define LOCAL_OFFSET atib + NCL_PER_SUPERCL * CL_SIZE
192 #ifndef REDUCE_SHUFFLE
193 /* shmem j force buffer */
194 __local float *f_buf = (__local float *)(LOCAL_OFFSET);
196 #define LOCAL_OFFSET f_buf + CL_SIZE * CL_SIZE * 3
198 /* Local buffer used to implement __any warp vote function from CUDA.
199 volatile is used to avoid compiler optimizations for AMD builds. */
200 volatile __local uint *warp_any = (__local uint*)(LOCAL_OFFSET);
203 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
204 sci = nb_sci.sci; /* super-cluster */
205 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
206 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
208 /* Pre-load i-atom x and q into shared memory */
209 ci = sci * NCL_PER_SUPERCL + tidxj;
210 ai = ci * CL_SIZE + tidxi;
212 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);
213 xqbuf.w *= nbparam->epsfac;
214 xqib[tidxj * CL_SIZE + tidxi] = xqbuf;
216 #ifdef IATYPE_SHMEM //NOTE: Should not be defined. Re-evaluate the effect of preloading at a suitable time.
217 /* Pre-load the i-atom types into shared memory */
218 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
220 /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
221 if(tidx==0 || tidx==32)
224 barrier(CLK_LOCAL_MEM_FENCE);
226 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
228 fci_buf[ci_offset] = (float3)(0.0f);
232 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
233 lje_coeff2 = nbparam->ewaldcoeff_lj*nbparam->ewaldcoeff_lj;
234 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
235 #endif /* LJ_EWALD */
242 #if defined EXCLUSION_FORCES /* Ewald or RF */
243 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
245 /* we have the diagonal: add the charge and LJ self interaction energy term */
246 for (i = 0; i < NCL_PER_SUPERCL; i++)
248 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
249 qi = xqib[i * CL_SIZE + tidxi].w;
253 E_lj += nbfp_climg2d[atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2];
254 #endif /* LJ_EWALD */
257 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
260 E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
261 #endif /* LJ_EWALD */
263 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
264 /* Correct for epsfac^2 due to adding qi^2 */
265 E_el /= nbparam->epsfac*CL_SIZE;
266 #if defined EL_RF || defined EL_CUTOFF
269 E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
271 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
273 #endif /* EXCLUSION_FORCES */
275 #endif /* CALC_ENERGIES */
277 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
278 for (j4 = cij4_start; j4 < cij4_end; j4++)
280 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
281 imask = pl_cj4[j4].imei[widx].imask;
282 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
288 /* Pre-load cj into shared memory on both warps separately */
289 if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
291 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
294 /* Unrolling this loop improves performance without pruning but
295 * with pruning it leads to slowdown.
297 * Tested with driver 1800.5
299 #if !defined PRUNE_NBL
303 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
305 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
307 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
309 cj = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
310 aj = cj * CL_SIZE + tidxj;
312 /* load j atom data */
314 xj = (float3)(xqbuf.xyz);
316 typej = atom_types[aj];
318 fcj_buf = (float3)(0.0f);
320 #if !defined PRUNE_NBL
323 for (i = 0; i < NCL_PER_SUPERCL; i++)
327 ci_offset = i; /* i force buffer offset */
329 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
330 ai = ci * CL_SIZE + tidxi; /* i atom index */
332 /* all threads load an atom from i cluster ci into shmem! */
333 xqbuf = xqib[i * CL_SIZE + tidxi];
334 xi = (float3)(xqbuf.xyz);
336 /* distance between i and j atoms */
341 /* vote.. should code shmem serialisation, wonder what the hit will be */
345 /* If _none_ of the atoms pairs are in cutoff range,
346 the bit corresponding to the current
347 cluster-pair in imask gets set to 0. */
355 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
357 /* cutoff & exclusion check */
358 #ifdef EXCLUSION_FORCES
359 if (r2 < rcoulomb_sq *
360 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
362 if (r2 < rcoulomb_sq * int_bit)
365 /* load the rest of the i-atom parameters */
367 #ifdef IATYPE_SHMEM //Should not be defined!
368 typei = atib[i * CL_SIZE + tidxi];
370 typei = atom_types[ai];
372 /* LJ 6*C6 and 12*C12 */
373 c6 = nbfp_climg2d[2 * (ntypes * typei + typej)];
374 c12 = nbfp_climg2d[2 * (ntypes * typei + typej)+1];
376 /* avoid NaN for excluded pairs at r=0 */
377 r2 += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
380 inv_r2 = inv_r * inv_r;
381 inv_r6 = inv_r2 * inv_r2 * inv_r2;
382 #if defined EXCLUSION_FORCES
383 /* We could mask inv_r2, but with Ewald
384 * masking both inv_r6 and F_invr is faster */
386 #endif /* EXCLUSION_FORCES */
388 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
389 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
390 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot)*ONE_TWELVETH_F -
391 c6 * (inv_r6 + nbparam->dispersion_shift.cpot)*ONE_SIXTH_F);
396 #ifdef LJ_FORCE_SWITCH
398 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
400 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
401 #endif /* CALC_ENERGIES */
402 #endif /* LJ_FORCE_SWITCH */
406 #ifdef LJ_EWALD_COMB_GEOM
408 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);
410 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
411 #endif /* CALC_ENERGIES */
412 #elif defined LJ_EWALD_COMB_LB
413 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
415 int_bit, true, &F_invr, &E_lj_p
418 #endif /* CALC_ENERGIES */
420 #endif /* LJ_EWALD_COMB_GEOM */
421 #endif /* LJ_EWALD */
423 #ifdef VDW_CUTOFF_CHECK
424 /* Separate VDW cut-off check to enable twin-range cut-offs
425 * (rvdw < rcoulomb <= rlist)
427 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
428 F_invr *= vdw_in_range;
430 E_lj_p *= vdw_in_range;
432 #endif /* VDW_CUTOFF_CHECK */
436 calculate_potential_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
438 calculate_potential_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
439 #endif /* CALC_ENERGIES */
440 #endif /* LJ_POT_SWITCH */
449 #ifdef EXCLUSION_FORCES
450 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
452 F_invr += qi * qj_f * inv_r2 * inv_r;
456 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
458 #if defined EL_EWALD_ANA
459 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
460 #elif defined EL_EWALD_TAB
461 F_invr += qi * qj_f * (int_bit*inv_r2 -
463 interpolate_coulomb_force_r(nbparam->coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
465 interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r, coulomb_tab_scale)
466 #endif /* USE_TEXOBJ */
468 #endif /* EL_EWALD_ANA/TAB */
472 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
475 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
478 /* 1.0f - erff is faster than erfcf */
479 E_el += qi * qj_f * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
480 #endif /* EL_EWALD_ANY */
484 /* accumulate j forces in registers */
487 /* accumulate i forces in registers */
488 fci_buf[ci_offset] += f_ij;
492 /* shift the mask bit by 1 */
496 /* reduce j forces */
498 /* store j forces in shmem */
499 f_buf[ tidx] = fcj_buf.x;
500 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
501 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
503 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
507 /* Update the imask with the new one which does not contain the
508 out of range clusters anymore. */
510 pl_cj4[j4].imei[widx].imask = imask;
515 /* skip central shifts when summing shift forces */
516 if (nb_sci.shift == CENTRAL)
523 /* reduce i forces */
524 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
526 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
528 f_buf[ tidx] = fci_buf[ci_offset].x;
529 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
530 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
531 barrier(CLK_LOCAL_MEM_FENCE);
532 reduce_force_i(f_buf, f,
533 &fshift_buf, bCalcFshift,
535 barrier(CLK_LOCAL_MEM_FENCE);
538 /* add up local shift forces into global mem */
541 /* Only threads with tidxj < 3 will update fshift.
542 The threads performing the update must be the same with the threads
543 which stored the reduction result in reduce_force_i function
546 atomicAdd_g_f(&(fshift[3 * nb_sci.shift + tidxj]), fshift_buf);
550 /* flush the energies to shmem and reduce them */
552 f_buf[FBUF_STRIDE + tidx] = E_el;
553 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
559 #undef EXCLUSION_FORCES