2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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 * \author Berk Hess <hess@kth.se>
45 * \ingroup module_mdlib
48 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
49 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/pbcutil/ishift.h"
52 /* Note that floating-point constants in CUDA code should be suffixed
53 * with f (e.g. 0.5f), to stop the compiler producing intermediate
54 * code that is in double precision.
57 #if GMX_PTX_ARCH < 300 && GMX_PTX_ARCH != 0
58 #error "nbnxn_cuda_kernel.cuh included with GMX_PTX_ARCH < 300 or host pass"
61 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
62 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
66 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
67 /* Macro to control the calculation of exclusion forces in the kernel
68 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
71 * Note: convenience macro, needs to be undef-ed at the end of the file.
73 #define EXCLUSION_FORCES
76 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
77 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
81 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
86 Kernel launch parameters:
87 - #blocks = #pair lists, blockId = pair list Id
88 - #threads = NTHREAD_Z * c_clSize^2
89 - shmem = see nbnxn_cuda.cu:calc_shmem_required_nonbonded()
91 Each thread calculates an i force-component taking one pair of i-j atoms.
95 /*! \brief Compute capability dependent definition of kernel launch configuration parameters.
97 * NTHREAD_Z controls the number of j-clusters processed concurrently on NTHREAD_Z
98 * warp-pairs per block.
100 * - On CC 2.0-3.5, and >=5.0 NTHREAD_Z == 1, translating to 64 th/block with 16
101 * blocks/multiproc, is the fastest even though this setup gives low occupancy
103 * NTHREAD_Z > 1 results in excessive register spilling unless the minimum blocks
104 * per multiprocessor is reduced proportionally to get the original number of max
105 * threads in flight (and slightly lower performance).
106 * - On CC 3.7 there are enough registers to double the number of threads; using
107 * NTHREADS_Z == 2 is fastest with 16 blocks (TODO: test with RF and other kernels
108 * with low-register use).
110 * Note that the current kernel implementation only supports NTHREAD_Z > 1 with
111 * shuffle-based reduction, hence CC >= 3.0.
114 * NOTEs on Volta / CUDA 9 extensions:
116 * - While active thread masks are required for the warp collectives
117 * (we use any and shfl), the kernel is designed such that all conditions
118 * (other than the inner-most distance check) including loop trip counts
119 * are warp-synchronous. Therefore, we don't need ballot to compute the
120 * active masks as these are all full-warp masks.
122 * - TODO: reconsider the use of __syncwarp(): its only role is currently to prevent
123 * WAR hazard due to the cj preload; we should try to replace it with direct
124 * loads (which may be faster given the improved L1 on Volta).
127 /* Kernel launch bounds for different compute capabilities. The value of NTHREAD_Z
128 * determines the number of threads per block and it is chosen such that
129 * 16 blocks/multiprocessor can be kept in flight.
130 * - CC 3.0,3.5, and >=5.0: NTHREAD_Z=1, (64, 16) bounds
131 * - CC 3.7: NTHREAD_Z=2, (128, 16) bounds
133 * Note: convenience macros, need to be undef-ed at the end of the file.
135 #if GMX_PTX_ARCH == 370
136 #define NTHREAD_Z (2)
137 #define MIN_BLOCKS_PER_MP (16)
139 #define NTHREAD_Z (1)
140 #define MIN_BLOCKS_PER_MP (16)
141 #endif /* GMX_PTX_ARCH == 370 */
142 #define THREADS_PER_BLOCK (c_clSize*c_clSize*NTHREAD_Z)
144 #if GMX_PTX_ARCH >= 350
145 #if (GMX_PTX_ARCH <= 210) && (NTHREAD_Z > 1)
146 #error NTHREAD_Z > 1 will give incorrect results on CC 2.x
149 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
151 __launch_bounds__(THREADS_PER_BLOCK)
152 #endif /* GMX_PTX_ARCH >= 350 */
155 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
157 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
158 #endif /* CALC_ENERGIES */
161 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
163 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
164 #endif /* CALC_ENERGIES */
165 #endif /* PRUNE_NBL */
166 (const cu_atomdata_t atdat,
167 const cu_nbparam_t nbparam,
168 const cu_plist_t plist,
170 #ifdef FUNCTION_DECLARATION_ONLY
171 ; /* Only do function declaration, omit the function body. */
174 /* convenience variables */
175 const nbnxn_sci_t *pl_sci = plist.sci;
179 nbnxn_cj4_t *pl_cj4 = plist.cj4;
180 const nbnxn_excl_t *excl = plist.excl;
182 const int *atom_types = atdat.atom_types;
183 int ntypes = atdat.ntypes;
185 const float2 *lj_comb = atdat.lj_comb;
186 float2 ljcp_i, ljcp_j;
188 const float4 *xq = atdat.xq;
190 const float3 *shift_vec = atdat.shift_vec;
191 float rcoulomb_sq = nbparam.rcoulomb_sq;
192 #ifdef VDW_CUTOFF_CHECK
193 float rvdw_sq = nbparam.rvdw_sq;
197 float lje_coeff2, lje_coeff6_6;
200 float two_k_rf = nbparam.two_k_rf;
203 float beta2 = nbparam.ewald_beta*nbparam.ewald_beta;
204 float beta3 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
207 float rlist_sq = nbparam.rlistOuter_sq;
212 float beta = nbparam.ewald_beta;
213 float ewald_shift = nbparam.sh_ewald;
215 float c_rf = nbparam.c_rf;
216 #endif /* EL_EWALD_ANY */
217 float *e_lj = atdat.e_lj;
218 float *e_el = atdat.e_el;
219 #endif /* CALC_ENERGIES */
221 /* thread/block/warp id-s */
222 unsigned int tidxi = threadIdx.x;
223 unsigned int tidxj = threadIdx.y;
224 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
226 unsigned int tidxz = 0;
228 unsigned int tidxz = threadIdx.z;
230 unsigned int bidx = blockIdx.x;
231 unsigned int widx = tidx / warp_size; /* warp index */
235 cij4_start, cij4_end;
239 int i, jm, j4, wexcl_idx;
242 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
243 float inv_r6, c6, c12;
246 float sigma, epsilon;
253 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
256 unsigned int wexcl, imask, mask_ji;
258 float3 xi, xj, rv, f_ij, fcj_buf;
259 float3 fci_buf[c_numClPerSupercl]; /* i force buffer */
262 /*! i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
263 const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
265 /* shmem buffer for i x+q pre-loading */
266 extern __shared__ float4 xqib[];
268 /* shmem buffer for cj, for each warp separately */
269 int *cjs = ((int *)(xqib + c_numClPerSupercl * c_clSize)) + tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
270 int *cjs_end = ((int *)(xqib + c_numClPerSupercl * c_clSize)) + NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
272 /* shmem buffer for i atom-type pre-loading */
275 /* shmem buffer for i-atom LJ combination rule parameters */
276 float2 *ljcpib = (float2 *)cjs_end;
279 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
280 sci = nb_sci.sci; /* super-cluster */
281 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
282 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
286 /* Pre-load i-atom x and q into shared memory */
287 ci = sci * c_numClPerSupercl + tidxj;
288 ai = ci * c_clSize + tidxi;
290 float *shiftptr = (float *)&shift_vec[nb_sci.shift];
291 xqbuf = xq[ai] + make_float4(LDG(shiftptr), LDG(shiftptr + 1), LDG(shiftptr + 2), 0.0f);
292 xqbuf.w *= nbparam.epsfac;
293 xqib[tidxj * c_clSize + tidxi] = xqbuf;
296 /* Pre-load the i-atom types into shared memory */
297 atib[tidxj * c_clSize + tidxi] = atom_types[ai];
299 /* Pre-load the LJ combination parameters into shared memory */
300 ljcpib[tidxj * c_clSize + tidxi] = lj_comb[ai];
305 for (i = 0; i < c_numClPerSupercl; i++)
307 fci_buf[i] = make_float3(0.0f);
311 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
312 lje_coeff2 = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
313 lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*c_oneSixth;
321 #ifdef EXCLUSION_FORCES /* Ewald or RF */
322 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*c_numClPerSupercl)
324 /* we have the diagonal: add the charge and LJ self interaction energy term */
325 for (i = 0; i < c_numClPerSupercl; i++)
327 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
328 qi = xqib[i * c_clSize + tidxi].w;
333 #if DISABLE_CUDA_TEXTURES
334 E_lj += LDG(&nbparam.nbfp[atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2]);
336 E_lj += tex1Dfetch<float>(nbparam.nbfp_texobj, atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2);
341 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
343 E_lj /= c_clSize*NTHREAD_Z;
344 E_lj *= 0.5f*c_oneSixth*lje_coeff6_6;
347 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
348 /* Correct for epsfac^2 due to adding qi^2 */
349 E_el /= nbparam.epsfac*c_clSize*NTHREAD_Z;
350 #if defined EL_RF || defined EL_CUTOFF
353 E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
355 #endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
357 #endif /* EXCLUSION_FORCES */
359 #endif /* CALC_ENERGIES */
361 #ifdef EXCLUSION_FORCES
362 const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
365 int j4LoopStart = cij4_start + tidxz;
366 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
367 for (j4 = j4LoopStart; j4 < cij4_end; j4 += NTHREAD_Z)
369 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
370 imask = pl_cj4[j4].imei[widx].imask;
371 wexcl = excl[wexcl_idx].pair[(tidx) & (warp_size - 1)];
377 /* Pre-load cj into shared memory on both warps separately */
378 if ((tidxj == 0 | tidxj == 4) & (tidxi < c_nbnxnGpuJgroupSize))
380 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize/c_splitClSize] = pl_cj4[j4].cj[tidxi];
382 gmx_syncwarp(c_fullWarpMask);
384 /* Unrolling this loop
385 - with pruning leads to register spilling;
386 - on Kepler and later it is much slower;
387 Tested with up to nvcc 7.5 */
388 for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
390 const unsigned int jmSkipCondition = imask & (superClInteractionMask << (jm * c_numClPerSupercl));
393 mask_ji = (1U << (jm * c_numClPerSupercl));
395 cj = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize/c_splitClSize];
396 aj = cj * c_clSize + tidxj;
398 /* load j atom data */
400 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
403 typej = atom_types[aj];
405 ljcp_j = lj_comb[aj];
408 fcj_buf = make_float3(0.0f);
410 #if !defined PRUNE_NBL
413 for (i = 0; i < c_numClPerSupercl; i++)
415 const unsigned int iInnerSkipCondition = imask & mask_ji;
416 if (iInnerSkipCondition)
418 ci = sci * c_numClPerSupercl + i; /* i cluster index */
420 /* all threads load an atom from i cluster ci into shmem! */
421 xqbuf = xqib[i * c_clSize + tidxi];
422 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
424 /* distance between i and j atoms */
429 /* If _none_ of the atoms pairs are in cutoff range,
430 the bit corresponding to the current
431 cluster-pair in imask gets set to 0. */
432 if (!gmx_any_sync(c_fullWarpMask, r2 < rlist_sq))
438 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
440 /* cutoff & exclusion check */
441 #ifdef EXCLUSION_FORCES
442 if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
444 if ((r2 < rcoulomb_sq) * int_bit)
447 /* load the rest of the i-atom parameters */
451 /* LJ 6*C6 and 12*C12 */
452 typei = atib[i * c_clSize + tidxi];
453 fetch_nbfp_c6_c12(c6, c12, nbparam, ntypes * typei + typej);
455 ljcp_i = ljcpib[i * c_clSize + tidxi];
457 c6 = ljcp_i.x * ljcp_j.x;
458 c12 = ljcp_i.y * ljcp_j.y;
460 /* LJ 2^(1/6)*sigma and 12*epsilon */
461 sigma = ljcp_i.x + ljcp_j.x;
462 epsilon = ljcp_i.y * ljcp_j.y;
463 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
464 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
466 #endif /* LJ_COMB_GEOM */
469 // Ensure distance do not become so small that r^-12 overflows
470 r2 = max(r2, NBNXN_MIN_RSQ);
473 inv_r2 = inv_r * inv_r;
474 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
475 inv_r6 = inv_r2 * inv_r2 * inv_r2;
476 #ifdef EXCLUSION_FORCES
477 /* We could mask inv_r2, but with Ewald
478 * masking both inv_r6 and F_invr is faster */
480 #endif /* EXCLUSION_FORCES */
482 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
483 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
484 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*c_oneTwelveth -
485 c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*c_oneSixth);
487 #else /* !LJ_COMB_LB || CALC_ENERGIES */
488 float sig_r = sigma*inv_r;
489 float sig_r2 = sig_r*sig_r;
490 float sig_r6 = sig_r2*sig_r2*sig_r2;
491 #ifdef EXCLUSION_FORCES
493 #endif /* EXCLUSION_FORCES */
495 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
496 #endif /* !LJ_COMB_LB || CALC_ENERGIES */
498 #ifdef LJ_FORCE_SWITCH
500 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
502 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
503 #endif /* CALC_ENERGIES */
504 #endif /* LJ_FORCE_SWITCH */
508 #ifdef LJ_EWALD_COMB_GEOM
510 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);
512 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
513 #endif /* CALC_ENERGIES */
514 #elif defined LJ_EWALD_COMB_LB
515 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
517 int_bit, &F_invr, &E_lj_p
520 #endif /* CALC_ENERGIES */
522 #endif /* LJ_EWALD_COMB_GEOM */
523 #endif /* LJ_EWALD */
527 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
529 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
530 #endif /* CALC_ENERGIES */
531 #endif /* LJ_POT_SWITCH */
533 #ifdef VDW_CUTOFF_CHECK
534 /* Separate VDW cut-off check to enable twin-range cut-offs
535 * (rvdw < rcoulomb <= rlist)
537 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
538 F_invr *= vdw_in_range;
540 E_lj_p *= vdw_in_range;
542 #endif /* VDW_CUTOFF_CHECK */
550 #ifdef EXCLUSION_FORCES
551 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
553 F_invr += qi * qj_f * inv_r2 * inv_r;
557 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
559 #if defined EL_EWALD_ANA
560 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
561 #elif defined EL_EWALD_TAB
562 F_invr += qi * qj_f * (int_bit*inv_r2 -
563 interpolate_coulomb_force_r(nbparam, r2 * inv_r)) * inv_r;
564 #endif /* EL_EWALD_ANA/TAB */
568 E_el += qi * qj_f * (int_bit*inv_r - c_rf);
571 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
574 /* 1.0f - erff is faster than erfcf */
575 E_el += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
576 #endif /* EL_EWALD_ANY */
580 /* accumulate j forces in registers */
583 /* accumulate i forces in registers */
588 /* shift the mask bit by 1 */
592 /* reduce j forces */
593 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj, c_fullWarpMask);
597 /* Update the imask with the new one which does not contain the
598 out of range clusters anymore. */
599 pl_cj4[j4].imei[widx].imask = imask;
602 // avoid shared memory WAR hazards between loop iterations
603 gmx_syncwarp(c_fullWarpMask);
606 /* skip central shifts when summing shift forces */
607 if (nb_sci.shift == CENTRAL)
612 float fshift_buf = 0.0f;
614 /* reduce i forces */
615 for (i = 0; i < c_numClPerSupercl; i++)
617 ai = (sci * c_numClPerSupercl + i) * c_clSize + tidxi;
618 reduce_force_i_warp_shfl(fci_buf[i], f,
619 &fshift_buf, bCalcFshift,
620 tidxj, ai, c_fullWarpMask);
623 /* add up local shift forces into global mem, tidxj indexes x,y,z */
624 if (bCalcFshift && (tidxj & 3) < 3)
626 atomicAdd(&(atdat.fshift[nb_sci.shift].x) + (tidxj & 3), fshift_buf);
630 /* reduce the energies over warps and store into global memory */
631 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx, c_fullWarpMask);
634 #endif /* FUNCTION_DECLARATION_ONLY */
637 #undef MIN_BLOCKS_PER_MP
638 #undef THREADS_PER_BLOCK
641 #undef EXCLUSION_FORCES