1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
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 */
50 Kernel launch parameters:
51 - #blocks = #pair lists, blockId = pair list Id
52 - #threads = CL_SIZE^2
53 - shmem = CL_SIZE^2 * sizeof(float)
55 Each thread calculates an i force-component taking one pair of i-j atoms.
59 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_prune)
61 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune)
65 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener)
67 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn)
70 (const cu_atomdata_t atdat,
71 const cu_nbparam_t nbparam,
72 const cu_plist_t plist,
75 /* convenience variables */
76 const nbnxn_sci_t *pl_sci = plist.sci;
80 nbnxn_cj4_t *pl_cj4 = plist.cj4;
81 const nbnxn_excl_t *excl = plist.excl;
82 const int *atom_types = atdat.atom_types;
83 int ntypes = atdat.ntypes;
84 const float4 *xq = atdat.xq;
86 const float3 *shift_vec = atdat.shift_vec;
87 float rcoulomb_sq = nbparam.rcoulomb_sq;
88 #ifdef VDW_CUTOFF_CHECK
89 float rvdw_sq = nbparam.rvdw_sq;
93 float two_k_rf = nbparam.two_k_rf;
96 float coulomb_tab_scale = nbparam.coulomb_tab_scale;
99 float rlist_sq = nbparam.rlist_sq;
103 float lj_shift = nbparam.sh_invrc6;
105 float beta = nbparam.ewald_beta;
106 float ewald_shift = nbparam.sh_ewald;
108 float c_rf = nbparam.c_rf;
110 float *e_lj = atdat.e_lj;
111 float *e_el = atdat.e_el;
114 /* thread/block/warp id-s */
115 unsigned int tidxi = threadIdx.x;
116 unsigned int tidxj = threadIdx.y;
117 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
118 unsigned int bidx = blockIdx.x;
119 unsigned int widx = tidx / WARP_SIZE; /* warp index */
121 int sci, ci, cj, ci_offset,
123 cij4_start, cij4_end,
125 i, jm, j4, wexcl_idx;
127 r2, inv_r, inv_r2, inv_r6,
134 unsigned int wexcl, imask, mask_ji;
136 float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
137 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
140 /* shmem buffer for i x+q pre-loading */
141 extern __shared__ float4 xqib[];
143 /* shmem buffer for i atom-type pre-loading */
144 int *atib = (int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
147 #ifndef REDUCE_SHUFFLE
148 /* shmem j force buffer */
150 float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
152 float *f_buf = (float *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
156 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
157 sci = nb_sci.sci; /* super-cluster */
158 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
159 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
161 /* Store the i-atom x and q in shared memory */
162 /* Note: the thread indexing here is inverted with respect to the
163 inner-loop as this results in slightly higher performance */
164 ci = sci * NCL_PER_SUPERCL + tidxi;
165 ai = ci * CL_SIZE + tidxj;
166 xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
168 ci = sci * NCL_PER_SUPERCL + tidxj;
169 ai = ci * CL_SIZE + tidxi;
170 atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
174 for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
176 fci_buf[ci_offset] = make_float3(0.0f);
183 #if defined EL_EWALD || defined EL_RF
184 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
186 /* we have the diagonal: add the charge self interaction energy term */
187 for (i = 0; i < NCL_PER_SUPERCL; i++)
189 qi = xqib[i * CL_SIZE + tidxi].w;
192 /* divide the self term equally over the j-threads */
195 E_el *= -nbparam.epsfac*0.5f*c_rf;
197 E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
203 /* skip central shifts when summing shift forces */
204 if (nb_sci.shift == CENTRAL)
209 fshift_buf = make_float3(0.0f);
211 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
212 for (j4 = cij4_start; j4 < cij4_end; j4++)
214 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
215 imask = pl_cj4[j4].imei[widx].imask;
216 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
222 /* Unrolling this loop
223 - with pruning leads to register spilling;
224 - on Kepler is much slower;
225 - doesn't work on CUDA <v4.1
226 Tested with nvcc 3.2 - 5.0.7 */
227 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && CUDA_VERSION >= 4010
230 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
232 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
234 mask_ji = (1U << (jm * NCL_PER_SUPERCL));
236 cj = pl_cj4[j4].cj[jm];
237 aj = cj * CL_SIZE + tidxj;
239 /* load j atom data */
241 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
242 qj_f = nbparam.epsfac * xqbuf.w;
243 typej = atom_types[aj];
245 fcj_buf = make_float3(0.0f);
247 /* The PME and RF kernels don't unroll with CUDA <v4.1. */
248 #if !defined PRUNE_NBL && !(CUDA_VERSION < 4010 && (defined EL_EWALD || defined EL_RF))
251 for(i = 0; i < NCL_PER_SUPERCL; i++)
255 ci_offset = i; /* i force buffer offset */
257 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
258 ai = ci * CL_SIZE + tidxi; /* i atom index */
260 /* all threads load an atom from i cluster ci into shmem! */
261 xqbuf = xqib[i * CL_SIZE + tidxi];
262 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
264 /* distance between i and j atoms */
269 /* If _none_ of the atoms pairs are in cutoff range,
270 the bit corresponding to the current
271 cluster-pair in imask gets set to 0. */
272 if (!__any(r2 < rlist_sq))
278 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
280 /* cutoff & exclusion check */
281 #if defined EL_EWALD || defined EL_RF
282 if (r2 < rcoulomb_sq *
283 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
285 if (r2 < rcoulomb_sq * int_bit)
288 /* load the rest of the i-atom parameters */
291 typei = atib[i * CL_SIZE + tidxi];
293 typei = atom_types[ai];
296 /* LJ 6*C6 and 12*C12 */
297 c6 = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej));
298 c12 = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej) + 1);
300 /* avoid NaN for excluded pairs at r=0 */
301 r2 += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
304 inv_r2 = inv_r * inv_r;
305 inv_r6 = inv_r2 * inv_r2 * inv_r2;
306 #if defined EL_EWALD || defined EL_RF
307 /* We could mask inv_r2, but with Ewald
308 * masking both inv_r6 and F_invr is faster */
312 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
315 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
318 #ifdef VDW_CUTOFF_CHECK
319 /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
320 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
321 F_invr *= vdw_in_range;
323 E_lj_p *= vdw_in_range;
332 F_invr += qi * qj_f * inv_r2 * inv_r;
335 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
338 F_invr += qi * qj_f * (int_bit*inv_r2 - interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)) * inv_r;
343 E_el += qi * qj_f * (inv_r - c_rf);
346 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
349 /* 1.0f - erff is faster than erfcf */
350 E_el += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
355 /* accumulate j forces in registers */
358 /* accumulate i forces in registers */
359 fci_buf[ci_offset] += f_ij;
363 /* shift the mask bit by 1 */
367 /* reduce j forces */
368 #ifdef REDUCE_SHUFFLE
369 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
371 /* store j forces in shmem */
372 f_buf[ tidx] = fcj_buf.x;
373 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
374 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
376 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
381 /* Update the imask with the new one which does not contain the
382 out of range clusters anymore. */
383 pl_cj4[j4].imei[widx].imask = imask;
388 /* reduce i forces */
389 for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
391 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
392 #ifdef REDUCE_SHUFFLE
393 reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
394 &fshift_buf, bCalcFshift,
397 f_buf[ tidx] = fci_buf[ci_offset].x;
398 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
399 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
401 reduce_force_i(f_buf, f,
402 &fshift_buf, bCalcFshift,
408 /* add up local shift forces into global mem */
409 #ifdef REDUCE_SHUFFLE
410 if (bCalcFshift && (tidxj == 0 || tidxj == 4))
412 if (bCalcFshift && tidxj == 0)
415 atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
416 atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
417 atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
421 #ifdef REDUCE_SHUFFLE
422 /* reduce the energies over warps and store into global memory */
423 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
425 /* flush the energies to shmem and reduce them */
427 f_buf[FBUF_STRIDE + tidx] = E_el;
428 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);