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.
43 Kernel launch parameters:
44 - #blocks = #pair lists, blockId = pair list Id
45 - #threads = CL_SIZE^2
46 - shmem = CL_SIZE^2 * sizeof(float)
48 Each thread calculates an i force-component taking one pair of i-j atoms.
52 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_prune_legacy)
54 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune_legacy)
58 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_legacy)
60 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _legacy)
63 (const cu_atomdata_t atdat,
64 const cu_nbparam_t nbparam,
65 const cu_plist_t plist,
68 /* convenience variables */
69 const nbnxn_sci_t *pl_sci = plist.sci;
73 nbnxn_cj4_t *pl_cj4 = plist.cj4;
74 const nbnxn_excl_t *excl = plist.excl;
75 const int *atom_types = atdat.atom_types;
76 int ntypes = atdat.ntypes;
77 const float4 *xq = atdat.xq;
79 const float3 *shift_vec = atdat.shift_vec;
80 float rcoulomb_sq = nbparam.rcoulomb_sq;
81 #ifdef VDW_CUTOFF_CHECK
82 float rvdw_sq = nbparam.rvdw_sq;
86 float two_k_rf = nbparam.two_k_rf;
89 float coulomb_tab_scale = nbparam.coulomb_tab_scale;
92 float rlist_sq = nbparam.rlist_sq;
96 float lj_shift = nbparam.sh_invrc6;
98 float beta = nbparam.ewald_beta;
99 float ewald_shift = nbparam.sh_ewald;
101 float c_rf = nbparam.c_rf;
103 float *e_lj = atdat.e_lj;
104 float *e_el = atdat.e_el;
107 /* thread/block/warp id-s */
108 unsigned int tidxi = threadIdx.x;
109 unsigned int tidxj = threadIdx.y;
110 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
111 unsigned int bidx = blockIdx.x;
112 unsigned int widx = tidx / WARP_SIZE; /* warp index */
114 int sci, ci, cj, ci_offset,
116 cij4_start, cij4_end,
118 i, cii, jm, j4, nsubi, wexcl_idx;
120 r2, inv_r, inv_r2, inv_r6,
126 unsigned int wexcl, int_bit, imask, imask_j;
128 unsigned int imask_prune;
131 float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
132 float3 fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
135 /* shmem buffer for i x+q pre-loading */
136 extern __shared__ float4 xqib[];
137 /* shmem j force buffer */
138 float *f_buf = (float *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
140 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
141 sci = nb_sci.sci; /* super-cluster */
142 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
143 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
145 /* Store the i-atom x and q in shared memory */
146 /* Note: the thread indexing here is inverted with respect to the
147 inner-loop as this results in slightly higher performance */
148 ci = sci * NCL_PER_SUPERCL + tidxi;
149 ai = ci * CL_SIZE + tidxj;
150 xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
153 for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
155 fci_buf[ci_offset] = make_float3(0.0f);
162 #if defined EL_EWALD || defined EL_RF
163 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
165 /* we have the diagonal: add the charge self interaction energy term */
166 for (i = 0; i < NCL_PER_SUPERCL; i++)
168 qi = xqib[i * CL_SIZE + tidxi].w;
171 /* divide the self term equally over the j-threads */
174 E_el *= -nbparam.epsfac*0.5f*c_rf;
176 E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
182 /* skip central shifts when summing shift forces */
183 if (nb_sci.shift == CENTRAL)
188 fshift_buf = make_float3(0.0f);
190 /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
191 for (j4 = cij4_start; j4 < cij4_end; j4++)
193 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
194 imask = pl_cj4[j4].imei[widx].imask;
195 wexcl = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
205 /* nvcc >v4.1 doesn't like this loop, it refuses to unroll it */
206 #if CUDA_VERSION >= 4010
209 for (jm = 0; jm < 4; jm++)
211 imask_j = (imask >> (jm * 8)) & 255U;
214 nsubi = __popc(imask_j);
216 cj = pl_cj4[j4].cj[jm];
217 aj = cj * CL_SIZE + tidxj;
219 /* load j atom data */
221 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
222 qj_f = nbparam.epsfac * xqbuf.w;
223 typej = atom_types[aj];
225 fcj_buf = make_float3(0.0f);
227 /* loop over the i-clusters in sci */
229 -- nvcc doesn't like my code, it refuses to unroll it
230 which is a pity because here unrolling could help. */
231 for (cii = 0; cii < nsubi; cii++)
233 i = __ffs(imask_j) - 1;
234 imask_j &= ~(1U << i);
236 ci_offset = i; /* i force buffer offset */
238 ci = sci * NCL_PER_SUPERCL + i; /* i cluster index */
239 ai = ci * CL_SIZE + tidxi; /* i atom index */
241 /* all threads load an atom from i cluster ci into shmem! */
242 xqbuf = xqib[i * CL_SIZE + tidxi];
243 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
245 /* distance between i and j atoms */
250 /* If _none_ of the atoms pairs are in cutoff range,
251 the bit corresponding to the current
252 cluster-pair in imask gets set to 0. */
253 if (!__any(r2 < rlist_sq))
255 imask_prune &= ~(1U << (jm * NCL_PER_SUPERCL + i));
259 int_bit = ((wexcl >> (jm * NCL_PER_SUPERCL + i)) & 1);
261 /* cutoff & exclusion check */
262 #if defined EL_EWALD || defined EL_RF
263 if (r2 < rcoulomb_sq *
264 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
266 if (r2 < rcoulomb_sq * int_bit)
269 /* load the rest of the i-atom parameters */
271 typei = atom_types[ai];
273 /* LJ 6*C6 and 12*C12 */
274 c6 = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej));
275 c12 = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej) + 1);
277 /* avoid NaN for excluded pairs at r=0 */
278 r2 += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
281 inv_r2 = inv_r * inv_r;
282 inv_r6 = inv_r2 * inv_r2 * inv_r2;
283 #if defined EL_EWALD || defined EL_RF
284 /* We could mask inv_r2, but with Ewald
285 * masking both inv_r6 and F_invr is faster */
289 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
292 E_lj_p = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
295 #ifdef VDW_CUTOFF_CHECK
296 /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
297 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
298 F_invr *= vdw_in_range;
300 E_lj_p *= vdw_in_range;
309 F_invr += qi * qj_f * inv_r2 * inv_r;
312 F_invr += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
315 F_invr += qi * qj_f * (int_bit*inv_r2 - interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)) * inv_r;
320 E_el += qi * qj_f * (inv_r - c_rf);
323 E_el += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
326 /* 1.0f - erff is faster than erfcf */
327 E_el += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
332 /* accumulate j forces in registers */
335 /* accumulate i forces in registers */
336 fci_buf[ci_offset] += f_ij;
340 /* store j forces in shmem */
341 f_buf[ tidx] = fcj_buf.x;
342 f_buf[ FBUF_STRIDE + tidx] = fcj_buf.y;
343 f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
345 /* reduce j forces */
346 reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
350 /* Update the imask with the new one which does not contain the
351 out of range clusters anymore. */
352 pl_cj4[j4].imei[widx].imask = imask_prune;
357 /* reduce i forces */
358 for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
360 ai = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
361 f_buf[ tidx] = fci_buf[ci_offset].x;
362 f_buf[ FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
363 f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
365 reduce_force_i(f_buf, f,
366 &fshift_buf, bCalcFshift,
371 /* add up local shift forces into global mem */
372 if (bCalcFshift && tidxj == 0)
374 atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
375 atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
376 atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
380 /* flush the energies to shmem and reduce them */
382 f_buf[FBUF_STRIDE + tidx] = E_el;
383 reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);