2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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.
38 * CUDA non-bonded kernel used through preprocessor-based code generation
39 * of multiple kernel flavors, see nbnxn_cuda_kernels.cuh.
41 * NOTE: No include fence as it is meant to be included multiple times.
43 * \author Szilárd Páll <pall.szilard@gmail.com>
44 * \ingroup module_mdlib
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/pbcutil/ishift.h"
50 /* Note that floating-point constants in CUDA code should be suffixed
51 * with f (e.g. 0.5f), to stop the compiler producing intermediate
52 * code that is in double precision.
55 #if __CUDA_ARCH__ >= 300
56 #define REDUCE_SHUFFLE
57 /* On Kepler pre-loading i-atom types to shmem gives a few %,
58 but on Fermi it does not */
62 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
63 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
67 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
68 /* Macro to control the calculation of exclusion forces in the kernel
69 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
72 * Note: convenience macro, needs to be undef-ed at the end of the file.
74 #define EXCLUSION_FORCES
77 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
78 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
83 Kernel launch parameters:
84 - #blocks = #pair lists, blockId = pair list Id
85 - #threads = NTHREAD_Z * CL_SIZE^2
86 - shmem = see nbnxn_cuda.cu:calc_shmem_required()
88 Each thread calculates an i force-component taking one pair of i-j atoms.
91 #if __CUDA_ARCH__ >= 350
92 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
94 __launch_bounds__(THREADS_PER_BLOCK)
98 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
100 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
104 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
106 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
109 (const cu_atomdata_t atdat,
110 const cu_nbparam_t nbparam,
111 const cu_plist_t plist,
114 /* convenience variables */
115 const nbnxn_sci_t *pl_sci = plist.sci;
119 nbnxn_cj4_t *pl_cj4 = plist.cj4;
120 const nbnxn_excl_t *excl = plist.excl;
121 const int *atom_types = atdat.atom_types;
122 int ntypes = atdat.ntypes;
123 const float4 *xq = atdat.xq;
125 const float3 *shift_vec = atdat.shift_vec;
126 float rcoulomb_sq = nbparam.rcoulomb_sq;
127 #ifdef VDW_CUTOFF_CHECK
128 float rvdw_sq = nbparam.rvdw_sq;
132 float lje_coeff2, lje_coeff6_6;
135 float two_k_rf = nbparam.two_k_rf;
138 float coulomb_tab_scale = nbparam.coulomb_tab_scale;
141 float beta2 = nbparam.ewald_beta*nbparam.ewald_beta;
142 float beta3 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
145 float rlist_sq = nbparam.rlist_sq;
150 float beta = nbparam.ewald_beta;
151 float ewald_shift = nbparam.sh_ewald;
153 float c_rf = nbparam.c_rf;
154 #endif /* EL_EWALD_ANY */
155 float *e_lj = atdat.e_lj;
156 float *e_el = atdat.e_el;
157 #endif /* CALC_ENERGIES */
159 /* thread/block/warp id-s */
160 unsigned int tidxi = threadIdx.x;
161 unsigned int tidxj = threadIdx.y;
162 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
164 unsigned int tidxz = 0;
166 unsigned int tidxz = threadIdx.z;
168 unsigned int bidx = blockIdx.x;
169 unsigned int widx = tidx / WARP_SIZE; /* warp index */
171 int sci, ci, cj, ci_offset,
173 cij4_start, cij4_end,
175 i, jm, j4, wexcl_idx;
177 r2, inv_r, inv_r2, inv_r6,
184 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
187 unsigned int wexcl, imask, mask_ji;
189 float3 xi, xj, rv, f_ij, fcj_buf;
190 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
193 /* shmem buffer for i x+q pre-loading */
194 extern __shared__ float4 xqib[];
195 /* shmem buffer for cj, for each warp separately */
196 int *cjs = ((int *)(xqib + NCL_PER_SUPERCL * CL_SIZE)) + tidxz * 2 * NBNXN_GPU_JGROUP_SIZE;
198 /* shmem buffer for i atom-type pre-loading */
199 int *atib = ((int *)(xqib + NCL_PER_SUPERCL * CL_SIZE)) + NTHREAD_Z * 2 * NBNXN_GPU_JGROUP_SIZE;
202 #ifndef REDUCE_SHUFFLE
203 /* shmem j force buffer */
205 float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
207 float *f_buf = (float *)(cjs + NTHREAD_Z * 2 * NBNXN_GPU_JGROUP_SIZE);
211 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
212 sci = nb_sci.sci; /* super-cluster */
213 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
214 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
218 /* Pre-load i-atom x and q into shared memory */
219 ci = sci * NCL_PER_SUPERCL + tidxj;
220 ai = ci * CL_SIZE + tidxi;
221 xqib[tidxj * CL_SIZE + tidxi] = xq[ai] + shift_vec[nb_sci.shift];
223 /* Pre-load the i-atom types into shared memory */
224 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
229 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
231 fci_buf[ci_offset] = make_float3(0.0f);
235 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
236 lje_coeff2 = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
237 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
238 #endif /* LJ_EWALD */
245 #if defined EXCLUSION_FORCES /* Ewald or RF */
246 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
248 /* we have the diagonal: add the charge and LJ self interaction energy term */
249 for (i = 0; i < NCL_PER_SUPERCL; i++)
251 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
252 qi = xqib[i * CL_SIZE + tidxi].w;
258 E_lj += tex1Dfetch<float>(nbparam.nbfp_texobj, atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2);
260 E_lj += tex1Dfetch(nbfp_texref, atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2);
261 #endif /* USE_TEXOBJ */
262 #endif /* LJ_EWALD */
266 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
268 E_lj /= CL_SIZE*NTHREAD_Z;
269 E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
270 #endif /* LJ_EWALD */
272 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
273 E_el /= CL_SIZE*NTHREAD_Z;
274 #if defined EL_RF || defined EL_CUTOFF
275 E_el *= -nbparam.epsfac*0.5f*c_rf;
277 E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
279 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
281 #endif /* EXCLUSION_FORCES */
283 #endif /* CALC_ENERGIES */
285 /* skip central shifts when summing shift forces */
286 if (nb_sci.shift == CENTRAL)
291 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
292 for (j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
294 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
295 imask = pl_cj4[j4].imei[widx].imask;
296 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
302 /* Pre-load cj into shared memory on both warps separately */
303 if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
305 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
308 /* Unrolling this loop
309 - with pruning leads to register spilling;
310 - on Kepler is much slower;
311 - doesn't work on CUDA <v4.1
312 Tested with nvcc 3.2 - 5.0.7 */
313 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && GMX_CUDA_VERSION >= 4010
316 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
318 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
320 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
322 cj = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
323 aj = cj * CL_SIZE + tidxj;
325 /* load j atom data */
327 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
328 qj_f = nbparam.epsfac * xqbuf.w;
329 typej = atom_types[aj];
331 fcj_buf = make_float3(0.0f);
333 /* The PME and RF kernels don't unroll with CUDA <v4.1. */
334 #if !defined PRUNE_NBL && !(GMX_CUDA_VERSION < 4010 && defined EXCLUSION_FORCES)
337 for (i = 0; i < NCL_PER_SUPERCL; i++)
341 ci_offset = i; /* i force buffer offset */
343 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
344 ai = ci * CL_SIZE + tidxi; /* i atom index */
346 /* all threads load an atom from i cluster ci into shmem! */
347 xqbuf = xqib[i * CL_SIZE + tidxi];
348 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
350 /* distance between i and j atoms */
355 /* If _none_ of the atoms pairs are in cutoff range,
356 the bit corresponding to the current
357 cluster-pair in imask gets set to 0. */
358 if (!__any(r2 < rlist_sq))
364 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
366 /* cutoff & exclusion check */
367 #ifdef EXCLUSION_FORCES
368 if (r2 < rcoulomb_sq *
369 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
371 if (r2 < rcoulomb_sq * int_bit)
374 /* load the rest of the i-atom parameters */
377 typei = atib[i * CL_SIZE + tidxi];
379 typei = atom_types[ai];
382 /* LJ 6*C6 and 12*C12 */
384 c6 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej));
385 c12 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej) + 1);
387 c6 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
388 c12 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
389 #endif /* USE_TEXOBJ */
392 /* avoid NaN for excluded pairs at r=0 */
393 r2 += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
396 inv_r2 = inv_r * inv_r;
397 inv_r6 = inv_r2 * inv_r2 * inv_r2;
398 #if defined EXCLUSION_FORCES
399 /* We could mask inv_r2, but with Ewald
400 * masking both inv_r6 and F_invr is faster */
402 #endif /* EXCLUSION_FORCES */
404 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
405 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
406 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*ONE_TWELVETH_F -
407 c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*ONE_SIXTH_F);
410 #ifdef LJ_FORCE_SWITCH
412 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
414 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
415 #endif /* CALC_ENERGIES */
416 #endif /* LJ_FORCE_SWITCH */
420 #ifdef LJ_EWALD_COMB_GEOM
422 calculate_lj_ewald_comb_geom_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
424 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
425 #endif /* CALC_ENERGIES */
426 #elif defined LJ_EWALD_COMB_LB
427 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
429 int_bit, &F_invr, &E_lj_p
432 #endif /* CALC_ENERGIES */
434 #endif /* LJ_EWALD_COMB_GEOM */
435 #endif /* LJ_EWALD */
437 #ifdef VDW_CUTOFF_CHECK
438 /* Separate VDW cut-off check to enable twin-range cut-offs
439 * (rvdw < rcoulomb <= rlist)
441 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
442 F_invr *= vdw_in_range;
444 E_lj_p *= vdw_in_range;
446 #endif /* VDW_CUTOFF_CHECK */
450 calculate_potential_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
452 calculate_potential_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
453 #endif /* CALC_ENERGIES */
454 #endif /* LJ_POT_SWITCH */
462 #ifdef EXCLUSION_FORCES
463 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
465 F_invr += qi * qj_f * inv_r2 * inv_r;
469 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
471 #if defined EL_EWALD_ANA
472 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
473 #elif defined EL_EWALD_TAB
474 F_invr += qi * qj_f * (int_bit*inv_r2 -
476 interpolate_coulomb_force_r(nbparam.coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
478 interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)
479 #endif /* USE_TEXOBJ */
481 #endif /* EL_EWALD_ANA/TAB */
485 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
488 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
491 /* 1.0f - erff is faster than erfcf */
492 E_el += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
493 #endif /* EL_EWALD_ANY */
497 /* accumulate j forces in registers */
500 /* accumulate i forces in registers */
501 fci_buf[ci_offset] += f_ij;
505 /* shift the mask bit by 1 */
509 /* reduce j forces */
510 #ifdef REDUCE_SHUFFLE
511 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
513 /* store j forces in shmem */
514 f_buf[ tidx] = fcj_buf.x;
515 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
516 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
518 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
523 /* Update the imask with the new one which does not contain the
524 out of range clusters anymore. */
525 pl_cj4[j4].imei[widx].imask = imask;
530 float fshift_buf = 0.0f;
532 /* reduce i forces */
533 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
535 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
536 #ifdef REDUCE_SHUFFLE
537 reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
538 &fshift_buf, bCalcFshift,
541 f_buf[ tidx] = fci_buf[ci_offset].x;
542 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
543 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
545 reduce_force_i(f_buf, f,
546 &fshift_buf, bCalcFshift,
552 /* add up local shift forces into global mem, tidxj indexes x,y,z */
553 #ifdef REDUCE_SHUFFLE
554 if (bCalcFshift && (tidxj & 3) < 3)
556 atomicAdd(&(atdat.fshift[nb_sci.shift].x) + (tidxj & ~4), fshift_buf);
559 if (bCalcFshift && tidxj < 3)
561 atomicAdd(&(atdat.fshift[nb_sci.shift].x) + tidxj, fshift_buf);
566 #ifdef REDUCE_SHUFFLE
567 /* reduce the energies over warps and store into global memory */
568 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
570 /* flush the energies to shmem and reduce them */
572 f_buf[FBUF_STRIDE + tidx] = E_el;
573 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
579 #undef EXCLUSION_FORCES