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(k_nbnxn, _ener_prune)
69 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune)
73 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener)
75 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn)
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;
115 float lj_shift = nbparam.sh_invrc6;
117 float beta = nbparam.ewald_beta;
118 float ewald_shift = nbparam.sh_ewald;
120 float c_rf = nbparam.c_rf;
122 float *e_lj = atdat.e_lj;
123 float *e_el = atdat.e_el;
126 /* thread/block/warp id-s */
127 unsigned int tidxi = threadIdx.x;
128 unsigned int tidxj = threadIdx.y;
129 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
130 unsigned int bidx = blockIdx.x;
131 unsigned int widx = tidx / WARP_SIZE; /* warp index */
133 int sci, ci, cj, ci_offset,
135 cij4_start, cij4_end,
137 i, jm, j4, wexcl_idx;
139 r2, inv_r, inv_r2, inv_r6,
144 float E_lj, E_el, E_lj_p;
146 unsigned int wexcl, imask, mask_ji;
148 float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
149 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
152 /* shmem buffer for i x+q pre-loading */
153 extern __shared__ float4 xqib[];
154 /* shmem buffer for cj, for both warps separately */
155 int *cjs = (int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
157 /* shmem buffer for i atom-type pre-loading */
158 int *atib = (int *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
161 #ifndef REDUCE_SHUFFLE
162 /* shmem j force buffer */
164 float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
166 float *f_buf = (float *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
170 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
171 sci = nb_sci.sci; /* super-cluster */
172 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
173 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
175 /* Store the i-atom x and q in shared memory */
176 /* Note: the thread indexing here is inverted with respect to the
177 inner-loop as this results in slightly higher performance */
178 ci = sci * NCL_PER_SUPERCL + tidxi;
179 ai = ci * CL_SIZE + tidxj;
180 xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
182 ci = sci * NCL_PER_SUPERCL + tidxj;
183 ai = ci * CL_SIZE + tidxi;
184 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
188 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
190 fci_buf[ci_offset] = make_float3(0.0f);
197 #if defined EL_EWALD_ANY || defined EL_RF
198 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
200 /* we have the diagonal: add the charge self interaction energy term */
201 for (i = 0; i < NCL_PER_SUPERCL; i++)
203 qi = xqib[i * CL_SIZE + tidxi].w;
206 /* divide the self term equally over the j-threads */
209 E_el *= -nbparam.epsfac*0.5f*c_rf;
211 E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
217 /* skip central shifts when summing shift forces */
218 if (nb_sci.shift == CENTRAL)
223 fshift_buf = make_float3(0.0f);
225 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
226 for (j4 = cij4_start; j4 < cij4_end; j4++)
228 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
229 imask = pl_cj4[j4].imei[widx].imask;
230 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
236 /* Pre-load cj into shared memory on both warps separately */
237 if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
239 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
242 /* Unrolling this loop
243 - with pruning leads to register spilling;
244 - on Kepler is much slower;
245 - doesn't work on CUDA <v4.1
246 Tested with nvcc 3.2 - 5.0.7 */
247 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && CUDA_VERSION >= 4010
250 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
252 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
254 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
256 cj = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
257 aj = cj * CL_SIZE + tidxj;
259 /* load j atom data */
261 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
262 qj_f = nbparam.epsfac * xqbuf.w;
263 typej = atom_types[aj];
265 fcj_buf = make_float3(0.0f);
267 /* The PME and RF kernels don't unroll with CUDA <v4.1. */
268 #if !defined PRUNE_NBL && !(CUDA_VERSION < 4010 && (defined EL_EWALD_ANY || defined EL_RF))
271 for (i = 0; i < NCL_PER_SUPERCL; i++)
275 ci_offset = i; /* i force buffer offset */
277 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
278 ai = ci * CL_SIZE + tidxi; /* i atom index */
280 /* all threads load an atom from i cluster ci into shmem! */
281 xqbuf = xqib[i * CL_SIZE + tidxi];
282 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
284 /* distance between i and j atoms */
289 /* If _none_ of the atoms pairs are in cutoff range,
290 the bit corresponding to the current
291 cluster-pair in imask gets set to 0. */
292 if (!__any(r2 < rlist_sq))
298 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
300 /* cutoff & exclusion check */
301 #if defined EL_EWALD_ANY || defined EL_RF
302 if (r2 < rcoulomb_sq *
303 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
305 if (r2 < rcoulomb_sq * int_bit)
308 /* load the rest of the i-atom parameters */
311 typei = atib[i * CL_SIZE + tidxi];
313 typei = atom_types[ai];
316 /* LJ 6*C6 and 12*C12 */
318 c6 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej));
319 c12 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej) + 1);
321 c6 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
322 c12 = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
323 #endif /* USE_TEXOBJ */
326 /* avoid NaN for excluded pairs at r=0 */
327 r2 += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
330 inv_r2 = inv_r * inv_r;
331 inv_r6 = inv_r2 * inv_r2 * inv_r2;
332 #if defined EL_EWALD_ANY || defined EL_RF
333 /* We could mask inv_r2, but with Ewald
334 * masking both inv_r6 and F_invr is faster */
338 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
341 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
344 #ifdef VDW_CUTOFF_CHECK
345 /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
346 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
347 F_invr *= vdw_in_range;
349 E_lj_p *= vdw_in_range;
358 F_invr += qi * qj_f * inv_r2 * inv_r;
361 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
363 #if defined EL_EWALD_ANA
364 F_invr += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
365 #elif defined EL_EWALD_TAB
366 F_invr += qi * qj_f * (int_bit*inv_r2 -
368 interpolate_coulomb_force_r(nbparam.coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
370 interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)
371 #endif /* USE_TEXOBJ */
373 #endif /* EL_EWALD_ANA/TAB */
377 E_el += qi * qj_f * (inv_r - c_rf);
380 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
383 /* 1.0f - erff is faster than erfcf */
384 E_el += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
385 #endif /* EL_EWALD_ANY */
389 /* accumulate j forces in registers */
392 /* accumulate i forces in registers */
393 fci_buf[ci_offset] += f_ij;
397 /* shift the mask bit by 1 */
401 /* reduce j forces */
402 #ifdef REDUCE_SHUFFLE
403 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
405 /* store j forces in shmem */
406 f_buf[ tidx] = fcj_buf.x;
407 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
408 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
410 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
415 /* Update the imask with the new one which does not contain the
416 out of range clusters anymore. */
417 pl_cj4[j4].imei[widx].imask = imask;
422 /* reduce i forces */
423 for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
425 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
426 #ifdef REDUCE_SHUFFLE
427 reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
428 &fshift_buf, bCalcFshift,
431 f_buf[ tidx] = fci_buf[ci_offset].x;
432 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
433 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
435 reduce_force_i(f_buf, f,
436 &fshift_buf, bCalcFshift,
442 /* add up local shift forces into global mem */
443 #ifdef REDUCE_SHUFFLE
444 if (bCalcFshift && (tidxj == 0 || tidxj == 4))
446 if (bCalcFshift && tidxj == 0)
449 atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
450 atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
451 atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
455 #ifdef REDUCE_SHUFFLE
456 /* reduce the energies over warps and store into global memory */
457 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
459 /* flush the energies to shmem and reduce them */
461 f_buf[FBUF_STRIDE + tidx] = E_el;
462 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);