2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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 "gromacs/math/utilities.h"
37 /* Note that floating-point constants in CUDA code should be suffixed
38 * with f (e.g. 0.5f), to stop the compiler producing intermediate
39 * code that is in double precision.
42 #if __CUDA_ARCH__ >= 300
43 #define REDUCE_SHUFFLE
44 /* On Kepler pre-loading i-atom types to shmem gives a few %,
45 but on Fermi it does not */
49 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
50 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
55 Kernel launch parameters:
56 - #blocks = #pair lists, blockId = pair list Id
57 - #threads = CL_SIZE^2
58 - shmem = CL_SIZE^2 * sizeof(float)
60 Each thread calculates an i force-component taking one pair of i-j atoms.
62 #if __CUDA_ARCH__ >= 350
63 __launch_bounds__(64, 16)
67 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
69 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
73 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
75 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
78 (const cu_atomdata_t atdat,
79 const cu_nbparam_t nbparam,
80 const cu_plist_t plist,
83 /* convenience variables */
84 const nbnxn_sci_t *pl_sci = plist.sci;
88 nbnxn_cj4_t *pl_cj4 = plist.cj4;
89 const nbnxn_excl_t *excl = plist.excl;
90 const int *atom_types = atdat.atom_types;
91 int ntypes = atdat.ntypes;
92 const float4 *xq = atdat.xq;
94 const float3 *shift_vec = atdat.shift_vec;
95 float rcoulomb_sq = nbparam.rcoulomb_sq;
96 #ifdef VDW_CUTOFF_CHECK
97 float rvdw_sq = nbparam.rvdw_sq;
101 float two_k_rf = nbparam.two_k_rf;
104 float coulomb_tab_scale = nbparam.coulomb_tab_scale;
107 float beta2 = nbparam.ewald_beta*nbparam.ewald_beta;
108 float beta3 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
111 float rlist_sq = nbparam.rlist_sq;
116 float beta = nbparam.ewald_beta;
117 float ewald_shift = nbparam.sh_ewald;
119 float c_rf = nbparam.c_rf;
120 #endif /* EL_EWALD_ANY */
121 float *e_lj = atdat.e_lj;
122 float *e_el = atdat.e_el;
123 #endif /* CALC_ENERGIES */
125 /* thread/block/warp id-s */
126 unsigned int tidxi = threadIdx.x;
127 unsigned int tidxj = threadIdx.y;
128 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
129 unsigned int bidx = blockIdx.x;
130 unsigned int widx = tidx / WARP_SIZE; /* warp index */
132 int sci, ci, cj, ci_offset,
134 cij4_start, cij4_end,
136 i, jm, j4, wexcl_idx;
138 r2, inv_r, inv_r2, inv_r6,
145 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
148 unsigned int wexcl, imask, mask_ji;
150 float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
151 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
154 /* shmem buffer for i x+q pre-loading */
155 extern __shared__ float4 xqib[];
156 /* shmem buffer for cj, for both warps separately */
157 int *cjs = (int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
159 /* shmem buffer for i atom-type pre-loading */
160 int *atib = (int *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
163 #ifndef REDUCE_SHUFFLE
164 /* shmem j force buffer */
166 float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
168 float *f_buf = (float *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
172 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
173 sci = nb_sci.sci; /* super-cluster */
174 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
175 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
177 /* Store the i-atom x and q in shared memory */
178 /* Note: the thread indexing here is inverted with respect to the
179 inner-loop as this results in slightly higher performance */
180 ci = sci * NCL_PER_SUPERCL + tidxi;
181 ai = ci * CL_SIZE + tidxj;
182 xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
184 ci = sci * NCL_PER_SUPERCL + tidxj;
185 ai = ci * CL_SIZE + tidxi;
186 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
190 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
192 fci_buf[ci_offset] = make_float3(0.0f);
199 #if defined EL_EWALD_ANY || defined EL_RF
200 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
202 /* we have the diagonal: add the charge self interaction energy term */
203 for (i = 0; i < NCL_PER_SUPERCL; i++)
205 qi = xqib[i * CL_SIZE + tidxi].w;
208 /* divide the self term equally over the j-threads */
211 E_el *= -nbparam.epsfac*0.5f*c_rf;
213 E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
219 /* skip central shifts when summing shift forces */
220 if (nb_sci.shift == CENTRAL)
225 fshift_buf = make_float3(0.0f);
227 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
228 for (j4 = cij4_start; j4 < cij4_end; j4++)
230 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
231 imask = pl_cj4[j4].imei[widx].imask;
232 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
238 /* Pre-load cj into shared memory on both warps separately */
239 if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
241 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
244 /* Unrolling this loop
245 - with pruning leads to register spilling;
246 - on Kepler is much slower;
247 - doesn't work on CUDA <v4.1
248 Tested with nvcc 3.2 - 5.0.7 */
249 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && CUDA_VERSION >= 4010
252 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
254 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
256 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
258 cj = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
259 aj = cj * CL_SIZE + tidxj;
261 /* load j atom data */
263 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
264 qj_f = nbparam.epsfac * xqbuf.w;
265 typej = atom_types[aj];
267 fcj_buf = make_float3(0.0f);
269 /* The PME and RF kernels don't unroll with CUDA <v4.1. */
270 #if !defined PRUNE_NBL && !(CUDA_VERSION < 4010 && (defined EL_EWALD_ANY || defined EL_RF))
273 for (i = 0; i < NCL_PER_SUPERCL; i++)
277 ci_offset = i; /* i force buffer offset */
279 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
280 ai = ci * CL_SIZE + tidxi; /* i atom index */
282 /* all threads load an atom from i cluster ci into shmem! */
283 xqbuf = xqib[i * CL_SIZE + tidxi];
284 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
286 /* distance between i and j atoms */
291 /* If _none_ of the atoms pairs are in cutoff range,
292 the bit corresponding to the current
293 cluster-pair in imask gets set to 0. */
294 if (!__any(r2 < rlist_sq))
300 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
302 /* cutoff & exclusion check */
303 #if defined EL_EWALD_ANY || defined EL_RF
304 if (r2 < rcoulomb_sq *
305 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
307 if (r2 < rcoulomb_sq * int_bit)
310 /* load the rest of the i-atom parameters */
313 typei = atib[i * CL_SIZE + tidxi];
315 typei = atom_types[ai];
318 /* LJ 6*C6 and 12*C12 */
320 c6 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej));
321 c12 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej) + 1);
323 c6 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
324 c12 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
325 #endif /* USE_TEXOBJ */
328 /* avoid NaN for excluded pairs at r=0 */
329 r2 += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
332 inv_r2 = inv_r * inv_r;
333 inv_r6 = inv_r2 * inv_r2 * inv_r2;
334 #if defined EL_EWALD_ANY || defined EL_RF
335 /* We could mask inv_r2, but with Ewald
336 * masking both inv_r6 and F_invr is faster */
340 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
341 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
342 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*ONE_TWELVETH_F -
343 c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*ONE_SIXTH_F);
346 #ifdef LJ_FORCE_SWITCH
348 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
350 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
351 #endif /* CALC_ENERGIES */
352 #endif /* LJ_FORCE_SWITCH */
354 #ifdef VDW_CUTOFF_CHECK
355 /* Separate VDW cut-off check to enable twin-range cut-offs
356 * (rvdw < rcoulomb <= rlist)
358 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
359 F_invr *= vdw_in_range;
361 E_lj_p *= vdw_in_range;
363 #endif /* VDW_CUTOFF_CHECK */
367 calculate_potential_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
369 calculate_potential_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
370 #endif /* CALC_ENERGIES */
371 #endif /* LJ_POT_SWITCH */
379 F_invr += qi * qj_f * inv_r2 * inv_r;
382 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
384 #if defined EL_EWALD_ANA
385 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
386 #elif defined EL_EWALD_TAB
387 F_invr += qi * qj_f * (int_bit*inv_r2 -
389 interpolate_coulomb_force_r(nbparam.coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
391 interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)
392 #endif /* USE_TEXOBJ */
394 #endif /* EL_EWALD_ANA/TAB */
398 E_el += qi * qj_f * (inv_r - c_rf);
401 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
404 /* 1.0f - erff is faster than erfcf */
405 E_el += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
406 #endif /* EL_EWALD_ANY */
410 /* accumulate j forces in registers */
413 /* accumulate i forces in registers */
414 fci_buf[ci_offset] += f_ij;
418 /* shift the mask bit by 1 */
422 /* reduce j forces */
423 #ifdef REDUCE_SHUFFLE
424 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
426 /* store j forces in shmem */
427 f_buf[ tidx] = fcj_buf.x;
428 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
429 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
431 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
436 /* Update the imask with the new one which does not contain the
437 out of range clusters anymore. */
438 pl_cj4[j4].imei[widx].imask = imask;
443 /* reduce i forces */
444 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
446 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
447 #ifdef REDUCE_SHUFFLE
448 reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
449 &fshift_buf, bCalcFshift,
452 f_buf[ tidx] = fci_buf[ci_offset].x;
453 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
454 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
456 reduce_force_i(f_buf, f,
457 &fshift_buf, bCalcFshift,
463 /* add up local shift forces into global mem */
464 #ifdef REDUCE_SHUFFLE
465 if (bCalcFshift && (tidxj == 0 || tidxj == 4))
467 if (bCalcFshift && tidxj == 0)
470 atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
471 atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
472 atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
476 #ifdef REDUCE_SHUFFLE
477 /* reduce the energies over warps and store into global memory */
478 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
480 /* flush the energies to shmem and reduce them */
482 f_buf[FBUF_STRIDE + tidx] = E_el;
483 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);