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 /* Note: convenience macros, need to be undef-ed at the end of the file. */
57 #define REDUCE_SHUFFLE
58 /* On Kepler pre-loading i-atom types to shmem gives a few %,
59 but on Fermi it does not */
63 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
64 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
68 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
69 /* Macro to control the calculation of exclusion forces in the kernel
70 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
73 * Note: convenience macro, needs to be undef-ed at the end of the file.
75 #define EXCLUSION_FORCES
78 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
79 /* 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 = NTHREAD_Z * CL_SIZE^2
88 - shmem = see nbnxn_cuda.cu:calc_shmem_required()
90 Each thread calculates an i force-component taking one pair of i-j atoms.
93 #if __CUDA_ARCH__ >= 350
94 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
96 __launch_bounds__(THREADS_PER_BLOCK)
100 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
102 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
106 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
108 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
111 (const cu_atomdata_t atdat,
112 const cu_nbparam_t nbparam,
113 const cu_plist_t plist,
116 /* convenience variables */
117 const nbnxn_sci_t *pl_sci = plist.sci;
121 nbnxn_cj4_t *pl_cj4 = plist.cj4;
122 const nbnxn_excl_t *excl = plist.excl;
123 const int *atom_types = atdat.atom_types;
124 int ntypes = atdat.ntypes;
125 const float4 *xq = atdat.xq;
127 const float3 *shift_vec = atdat.shift_vec;
128 float rcoulomb_sq = nbparam.rcoulomb_sq;
129 #ifdef VDW_CUTOFF_CHECK
130 float rvdw_sq = nbparam.rvdw_sq;
134 float lje_coeff2, lje_coeff6_6;
137 float two_k_rf = nbparam.two_k_rf;
140 float coulomb_tab_scale = nbparam.coulomb_tab_scale;
143 float beta2 = nbparam.ewald_beta*nbparam.ewald_beta;
144 float beta3 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
147 float rlist_sq = nbparam.rlist_sq;
152 float beta = nbparam.ewald_beta;
153 float ewald_shift = nbparam.sh_ewald;
155 float c_rf = nbparam.c_rf;
156 #endif /* EL_EWALD_ANY */
157 float *e_lj = atdat.e_lj;
158 float *e_el = atdat.e_el;
159 #endif /* CALC_ENERGIES */
161 /* thread/block/warp id-s */
162 unsigned int tidxi = threadIdx.x;
163 unsigned int tidxj = threadIdx.y;
164 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
166 unsigned int tidxz = 0;
168 unsigned int tidxz = threadIdx.z;
170 unsigned int bidx = blockIdx.x;
171 unsigned int widx = tidx / WARP_SIZE; /* warp index */
173 int sci, ci, cj, ci_offset,
175 cij4_start, cij4_end,
177 i, jm, j4, wexcl_idx;
179 r2, inv_r, inv_r2, inv_r6,
186 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
189 unsigned int wexcl, imask, mask_ji;
191 float3 xi, xj, rv, f_ij, fcj_buf;
192 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
195 /* shmem buffer for i x+q pre-loading */
196 extern __shared__ float4 xqib[];
197 /* shmem buffer for cj, for each warp separately */
198 int *cjs = ((int *)(xqib + NCL_PER_SUPERCL * CL_SIZE)) + tidxz * 2 * NBNXN_GPU_JGROUP_SIZE;
200 /* shmem buffer for i atom-type pre-loading */
201 int *atib = ((int *)(xqib + NCL_PER_SUPERCL * CL_SIZE)) + NTHREAD_Z * 2 * NBNXN_GPU_JGROUP_SIZE;
204 #ifndef REDUCE_SHUFFLE
205 /* shmem j force buffer */
207 float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
209 float *f_buf = (float *)(cjs + NTHREAD_Z * 2 * NBNXN_GPU_JGROUP_SIZE);
213 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
214 sci = nb_sci.sci; /* super-cluster */
215 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
216 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
220 /* Pre-load i-atom x and q into shared memory */
221 ci = sci * NCL_PER_SUPERCL + tidxj;
222 ai = ci * CL_SIZE + tidxi;
223 xqib[tidxj * CL_SIZE + tidxi] = xq[ai] + shift_vec[nb_sci.shift];
225 /* Pre-load the i-atom types into shared memory */
226 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
231 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
233 fci_buf[ci_offset] = make_float3(0.0f);
237 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
238 lje_coeff2 = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
239 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
240 #endif /* LJ_EWALD */
247 #if defined EXCLUSION_FORCES /* Ewald or RF */
248 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
250 /* we have the diagonal: add the charge and LJ self interaction energy term */
251 for (i = 0; i < NCL_PER_SUPERCL; i++)
253 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
254 qi = xqib[i * CL_SIZE + tidxi].w;
260 E_lj += tex1Dfetch<float>(nbparam.nbfp_texobj, atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2);
262 E_lj += tex1Dfetch(nbfp_texref, atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2);
263 #endif /* USE_TEXOBJ */
264 #endif /* LJ_EWALD */
268 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
270 E_lj /= CL_SIZE*NTHREAD_Z;
271 E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
272 #endif /* LJ_EWALD */
274 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
275 E_el /= CL_SIZE*NTHREAD_Z;
276 #if defined EL_RF || defined EL_CUTOFF
277 E_el *= -nbparam.epsfac*0.5f*c_rf;
279 E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
281 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
283 #endif /* EXCLUSION_FORCES */
285 #endif /* CALC_ENERGIES */
287 /* skip central shifts when summing shift forces */
288 if (nb_sci.shift == CENTRAL)
293 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
294 for (j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
296 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
297 imask = pl_cj4[j4].imei[widx].imask;
298 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
304 /* Pre-load cj into shared memory on both warps separately */
305 if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
307 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
310 /* Unrolling this loop
311 - with pruning leads to register spilling;
312 - on Kepler is much slower;
313 - doesn't work on CUDA <v4.1
314 Tested with nvcc 3.2 - 5.0.7 */
315 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && GMX_CUDA_VERSION >= 4010
318 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
320 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
322 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
324 cj = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
325 aj = cj * CL_SIZE + tidxj;
327 /* load j atom data */
329 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
330 qj_f = nbparam.epsfac * xqbuf.w;
331 typej = atom_types[aj];
333 fcj_buf = make_float3(0.0f);
335 /* The PME and RF kernels don't unroll with CUDA <v4.1. */
336 #if !defined PRUNE_NBL && !(GMX_CUDA_VERSION < 4010 && defined EXCLUSION_FORCES)
339 for (i = 0; i < NCL_PER_SUPERCL; i++)
343 ci_offset = i; /* i force buffer offset */
345 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
346 ai = ci * CL_SIZE + tidxi; /* i atom index */
348 /* all threads load an atom from i cluster ci into shmem! */
349 xqbuf = xqib[i * CL_SIZE + tidxi];
350 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
352 /* distance between i and j atoms */
357 /* If _none_ of the atoms pairs are in cutoff range,
358 the bit corresponding to the current
359 cluster-pair in imask gets set to 0. */
360 if (!__any(r2 < rlist_sq))
366 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
368 /* cutoff & exclusion check */
369 #ifdef EXCLUSION_FORCES
370 if (r2 < rcoulomb_sq *
371 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
373 if (r2 < rcoulomb_sq * int_bit)
376 /* load the rest of the i-atom parameters */
379 typei = atib[i * CL_SIZE + tidxi];
381 typei = atom_types[ai];
384 /* LJ 6*C6 and 12*C12 */
386 c6 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej));
387 c12 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej) + 1);
389 c6 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
390 c12 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
391 #endif /* USE_TEXOBJ */
394 /* avoid NaN for excluded pairs at r=0 */
395 r2 += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
398 inv_r2 = inv_r * inv_r;
399 inv_r6 = inv_r2 * inv_r2 * inv_r2;
400 #if defined EXCLUSION_FORCES
401 /* We could mask inv_r2, but with Ewald
402 * masking both inv_r6 and F_invr is faster */
404 #endif /* EXCLUSION_FORCES */
406 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
407 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
408 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*ONE_TWELVETH_F -
409 c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*ONE_SIXTH_F);
412 #ifdef LJ_FORCE_SWITCH
414 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
416 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
417 #endif /* CALC_ENERGIES */
418 #endif /* LJ_FORCE_SWITCH */
422 #ifdef LJ_EWALD_COMB_GEOM
424 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);
426 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
427 #endif /* CALC_ENERGIES */
428 #elif defined LJ_EWALD_COMB_LB
429 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
431 int_bit, &F_invr, &E_lj_p
434 #endif /* CALC_ENERGIES */
436 #endif /* LJ_EWALD_COMB_GEOM */
437 #endif /* LJ_EWALD */
439 #ifdef VDW_CUTOFF_CHECK
440 /* Separate VDW cut-off check to enable twin-range cut-offs
441 * (rvdw < rcoulomb <= rlist)
443 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
444 F_invr *= vdw_in_range;
446 E_lj_p *= vdw_in_range;
448 #endif /* VDW_CUTOFF_CHECK */
452 calculate_potential_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
454 calculate_potential_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
455 #endif /* CALC_ENERGIES */
456 #endif /* LJ_POT_SWITCH */
464 #ifdef EXCLUSION_FORCES
465 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
467 F_invr += qi * qj_f * inv_r2 * inv_r;
471 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
473 #if defined EL_EWALD_ANA
474 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
475 #elif defined EL_EWALD_TAB
476 F_invr += qi * qj_f * (int_bit*inv_r2 -
478 interpolate_coulomb_force_r(nbparam.coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
480 interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)
481 #endif /* USE_TEXOBJ */
483 #endif /* EL_EWALD_ANA/TAB */
487 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
490 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
493 /* 1.0f - erff is faster than erfcf */
494 E_el += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
495 #endif /* EL_EWALD_ANY */
499 /* accumulate j forces in registers */
502 /* accumulate i forces in registers */
503 fci_buf[ci_offset] += f_ij;
507 /* shift the mask bit by 1 */
511 /* reduce j forces */
512 #ifdef REDUCE_SHUFFLE
513 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
515 /* store j forces in shmem */
516 f_buf[ tidx] = fcj_buf.x;
517 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
518 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
520 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
525 /* Update the imask with the new one which does not contain the
526 out of range clusters anymore. */
527 pl_cj4[j4].imei[widx].imask = imask;
532 float fshift_buf = 0.0f;
534 /* reduce i forces */
535 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
537 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
538 #ifdef REDUCE_SHUFFLE
539 reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
540 &fshift_buf, bCalcFshift,
543 f_buf[ tidx] = fci_buf[ci_offset].x;
544 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
545 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
547 reduce_force_i(f_buf, f,
548 &fshift_buf, bCalcFshift,
554 /* add up local shift forces into global mem, tidxj indexes x,y,z */
555 #ifdef REDUCE_SHUFFLE
556 if (bCalcFshift && (tidxj & 3) < 3)
558 atomicAdd(&(atdat.fshift[nb_sci.shift].x) + (tidxj & ~4), fshift_buf);
561 if (bCalcFshift && tidxj < 3)
563 atomicAdd(&(atdat.fshift[nb_sci.shift].x) + tidxj, fshift_buf);
568 #ifdef REDUCE_SHUFFLE
569 /* reduce the energies over warps and store into global memory */
570 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
572 /* flush the energies to shmem and reduce them */
574 f_buf[FBUF_STRIDE + tidx] = E_el;
575 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
580 #undef REDUCE_SHUFFLE
584 #undef EXCLUSION_FORCES