2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
39 * CUDA non-bonded kernel used through preprocessor-based code generation
40 * of multiple kernel flavors, see nbnxn_cuda_kernels.cuh.
42 * NOTE: No include fence as it is meant to be included multiple times.
44 * \author Szilárd Páll <pall.szilard@gmail.com>
45 * \author Berk Hess <hess@kth.se>
46 * \ingroup module_nbnxm
49 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
50 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
51 #include "gromacs/gpu_utils/typecasts.cuh"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/pbcutil/ishift.h"
55 /* Note that floating-point constants in CUDA code should be suffixed
56 * with f (e.g. 0.5f), to stop the compiler producing intermediate
57 * code that is in double precision.
60 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
61 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
65 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
66 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
70 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD \
71 || (defined EL_CUTOFF && defined CALC_ENERGIES)
72 /* Macro to control the calculation of exclusion forces in the kernel
73 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
76 * Note: convenience macro, needs to be undef-ed at the end of the file.
78 # define EXCLUSION_FORCES
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 3.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
146 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
148 __launch_bounds__(THREADS_PER_BLOCK)
149 #endif /* GMX_PTX_ARCH >= 350 */
151 # ifdef CALC_ENERGIES
152 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
154 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
155 # endif /* CALC_ENERGIES */
157 # ifdef CALC_ENERGIES
158 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
160 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
161 # endif /* CALC_ENERGIES */
162 #endif /* PRUNE_NBL */
163 (const NBAtomDataGpu atdat, const NBParamGpu nbparam, const Nbnxm::gpu_plist plist, bool bCalcFshift)
164 #ifdef FUNCTION_DECLARATION_ONLY
165 ; /* Only do function declaration, omit the function body. */
168 /* convenience variables */
169 const nbnxn_sci_t* pl_sci = plist.sci;
173 nbnxn_cj4_t* pl_cj4 = plist.cj4;
174 const nbnxn_excl_t* excl = plist.excl;
176 const int* atom_types = atdat.atomTypes;
177 int ntypes = atdat.numTypes;
179 const float2* lj_comb = atdat.ljComb;
180 float2 ljcp_i, ljcp_j;
182 const float4* xq = atdat.xq;
183 float3* f = asFloat3(atdat.f);
184 const float3* shift_vec = asFloat3(atdat.shiftVec);
185 float rcoulomb_sq = nbparam.rcoulomb_sq;
186 # ifdef VDW_CUTOFF_CHECK
187 float rvdw_sq = nbparam.rvdw_sq;
191 float lje_coeff2, lje_coeff6_6;
194 float two_k_rf = nbparam.two_k_rf;
197 float beta2 = nbparam.ewald_beta * nbparam.ewald_beta;
198 float beta3 = nbparam.ewald_beta * nbparam.ewald_beta * nbparam.ewald_beta;
201 float rlist_sq = nbparam.rlistOuter_sq;
204 # ifdef CALC_ENERGIES
206 float beta = nbparam.ewald_beta;
207 float ewald_shift = nbparam.sh_ewald;
209 float reactionFieldShift = nbparam.c_rf;
210 # endif /* EL_EWALD_ANY */
211 float* e_lj = atdat.eLJ;
212 float* e_el = atdat.eElec;
213 # endif /* CALC_ENERGIES */
215 /* thread/block/warp id-s */
216 unsigned int tidxi = threadIdx.x;
217 unsigned int tidxj = threadIdx.y;
218 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
220 unsigned int tidxz = 0;
222 unsigned int tidxz = threadIdx.z;
224 unsigned int bidx = blockIdx.x;
225 unsigned int widx = tidx / warp_size; /* warp index */
227 int sci, ci, cj, ai, aj, cij4_start, cij4_end;
231 int i, jm, j4, wexcl_idx;
232 float qi, qj_f, r2, inv_r, inv_r2;
233 # if !defined LJ_COMB_LB || defined CALC_ENERGIES
234 float inv_r6, c6, c12;
237 float sigma, epsilon;
239 float int_bit, F_invr;
240 # ifdef CALC_ENERGIES
243 # if defined CALC_ENERGIES || defined LJ_POT_SWITCH
246 unsigned int wexcl, imask, mask_ji;
248 float3 xi, xj, rv, f_ij, fcj_buf;
249 float3 fci_buf[c_nbnxnGpuNumClusterPerSupercluster]; /* i force buffer */
252 /*! i-cluster interaction mask for a super-cluster with all c_nbnxnGpuNumClusterPerSupercluster=8 bits set */
253 const unsigned superClInteractionMask = ((1U << c_nbnxnGpuNumClusterPerSupercluster) - 1U);
255 /*********************************************************************
256 * Set up shared memory pointers.
257 * sm_nextSlotPtr should always be updated to point to the "next slot",
258 * that is past the last point where data has been stored.
260 extern __shared__ char sm_dynamicShmem[];
261 char* sm_nextSlotPtr = sm_dynamicShmem;
262 static_assert(sizeof(char) == 1,
263 "The shared memory offset calculation assumes that char is 1 byte");
265 /* shmem buffer for i x+q pre-loading */
266 float4* xqib = (float4*)sm_nextSlotPtr;
267 sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*xqib));
269 /* shmem buffer for cj, for each warp separately */
270 int* cjs = (int*)(sm_nextSlotPtr);
271 /* the cjs buffer's use expects a base pointer offset for pairs of warps in the j-concurrent execution */
272 cjs += tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
273 sm_nextSlotPtr += (NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(*cjs));
276 /* shmem buffer for i atom-type pre-loading */
277 int* atib = (int*)sm_nextSlotPtr;
278 sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*atib));
280 /* shmem buffer for i-atom LJ combination rule parameters */
281 float2* ljcpib = (float2*)sm_nextSlotPtr;
282 sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*ljcpib));
284 /*********************************************************************/
286 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
287 sci = nb_sci.sci; /* super-cluster */
288 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
289 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
293 /* Pre-load i-atom x and q into shared memory */
294 ci = sci * c_nbnxnGpuNumClusterPerSupercluster + tidxj;
295 ai = ci * c_clSize + tidxi;
297 float* shiftptr = (float*)&shift_vec[nb_sci.shift];
298 xqbuf = xq[ai] + make_float4(LDG(shiftptr), LDG(shiftptr + 1), LDG(shiftptr + 2), 0.0f);
299 xqbuf.w *= nbparam.epsfac;
300 xqib[tidxj * c_clSize + tidxi] = xqbuf;
303 /* Pre-load the i-atom types into shared memory */
304 atib[tidxj * c_clSize + tidxi] = atom_types[ai];
306 /* Pre-load the LJ combination parameters into shared memory */
307 ljcpib[tidxj * c_clSize + tidxi] = lj_comb[ai];
312 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
314 fci_buf[i] = make_float3(0.0f);
318 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
319 lje_coeff2 = nbparam.ewaldcoeff_lj * nbparam.ewaldcoeff_lj;
320 lje_coeff6_6 = lje_coeff2 * lje_coeff2 * lje_coeff2 * c_oneSixth;
324 # ifdef CALC_ENERGIES
328 # ifdef EXCLUSION_FORCES /* Ewald or RF */
329 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci * c_nbnxnGpuNumClusterPerSupercluster)
331 /* we have the diagonal: add the charge and LJ self interaction energy term */
332 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
334 # if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
335 qi = xqib[i * c_clSize + tidxi].w;
340 // load only the first 4 bytes of the parameter pair (equivalent with nbfp[idx].x)
341 E_lj += LDG((float*)&nbparam.nbfp[atom_types[(sci * c_nbnxnGpuNumClusterPerSupercluster + i) * c_clSize + tidxi]
346 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
348 E_lj /= c_clSize * NTHREAD_Z;
349 E_lj *= 0.5f * c_oneSixth * lje_coeff6_6;
352 # if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
353 /* Correct for epsfac^2 due to adding qi^2 */
354 E_el /= nbparam.epsfac * c_clSize * NTHREAD_Z;
355 # if defined EL_RF || defined EL_CUTOFF
356 E_el *= -0.5f * reactionFieldShift;
358 E_el *= -beta * M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
360 # endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
362 # endif /* EXCLUSION_FORCES */
364 # endif /* CALC_ENERGIES */
366 # ifdef EXCLUSION_FORCES
367 const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
370 /* loop over the j clusters = seen by any of the atoms in the current super-cluster;
371 * The loop stride NTHREAD_Z ensures that consecutive warps-pairs are assigned
372 * consecutive j4's entries.
374 for (j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
376 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
377 imask = pl_cj4[j4].imei[widx].imask;
378 wexcl = excl[wexcl_idx].pair[(tidx) & (warp_size - 1)];
384 /* Pre-load cj into shared memory on both warps separately */
385 if ((tidxj == 0 | tidxj == 4) & (tidxi < c_nbnxnGpuJgroupSize))
387 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize / c_splitClSize] = pl_cj4[j4].cj[tidxi];
389 __syncwarp(c_fullWarpMask);
391 /* Unrolling this loop
392 - with pruning leads to register spilling;
393 - on Kepler and later it is much slower;
394 Tested with up to nvcc 7.5 */
395 for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
397 if (imask & (superClInteractionMask << (jm * c_nbnxnGpuNumClusterPerSupercluster)))
399 mask_ji = (1U << (jm * c_nbnxnGpuNumClusterPerSupercluster));
401 cj = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize / c_splitClSize];
402 aj = cj * c_clSize + tidxj;
404 /* load j atom data */
406 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
409 typej = atom_types[aj];
411 ljcp_j = lj_comb[aj];
414 fcj_buf = make_float3(0.0f);
416 # if !defined PRUNE_NBL
419 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
423 ci = sci * c_nbnxnGpuNumClusterPerSupercluster + i; /* i cluster index */
425 /* all threads load an atom from i cluster ci into shmem! */
426 xqbuf = xqib[i * c_clSize + tidxi];
427 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
429 /* distance between i and j atoms */
434 /* If _none_ of the atoms pairs are in cutoff range,
435 the bit corresponding to the current
436 cluster-pair in imask gets set to 0. */
437 if (!__any_sync(c_fullWarpMask, r2 < rlist_sq))
443 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
445 /* cutoff & exclusion check */
446 # ifdef EXCLUSION_FORCES
447 if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
449 if ((r2 < rcoulomb_sq) * int_bit)
452 /* load the rest of the i-atom parameters */
456 /* LJ 6*C6 and 12*C12 */
457 typei = atib[i * c_clSize + tidxi];
458 fetch_nbfp_c6_c12(c6, c12, nbparam, ntypes * typei + typej);
460 ljcp_i = ljcpib[i * c_clSize + tidxi];
462 c6 = ljcp_i.x * ljcp_j.x;
463 c12 = ljcp_i.y * ljcp_j.y;
465 /* LJ 2^(1/6)*sigma and 12*epsilon */
466 sigma = ljcp_i.x + ljcp_j.x;
467 epsilon = ljcp_i.y * ljcp_j.y;
468 # if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
469 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
471 # endif /* LJ_COMB_GEOM */
472 # endif /* LJ_COMB */
474 // Ensure distance do not become so small that r^-12 overflows
475 r2 = max(r2, c_nbnxnMinDistanceSquared);
478 inv_r2 = inv_r * inv_r;
479 # if !defined LJ_COMB_LB || defined CALC_ENERGIES
480 inv_r6 = inv_r2 * inv_r2 * inv_r2;
481 # ifdef EXCLUSION_FORCES
482 /* We could mask inv_r2, but with Ewald
483 * masking both inv_r6 and F_invr is faster */
485 # endif /* EXCLUSION_FORCES */
487 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
488 # if defined CALC_ENERGIES || defined LJ_POT_SWITCH
490 * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot) * c_oneTwelveth
491 - c6 * (inv_r6 + nbparam.dispersion_shift.cpot) * c_oneSixth);
493 # else /* !LJ_COMB_LB || CALC_ENERGIES */
494 float sig_r = sigma * inv_r;
495 float sig_r2 = sig_r * sig_r;
496 float sig_r6 = sig_r2 * sig_r2 * sig_r2;
497 # ifdef EXCLUSION_FORCES
499 # endif /* EXCLUSION_FORCES */
501 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
502 # endif /* !LJ_COMB_LB || CALC_ENERGIES */
504 # ifdef LJ_FORCE_SWITCH
505 # ifdef CALC_ENERGIES
506 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
508 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
509 # endif /* CALC_ENERGIES */
510 # endif /* LJ_FORCE_SWITCH */
514 # ifdef LJ_EWALD_COMB_GEOM
515 # ifdef CALC_ENERGIES
516 calculate_lj_ewald_comb_geom_F_E(
517 nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
519 calculate_lj_ewald_comb_geom_F(
520 nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
521 # endif /* CALC_ENERGIES */
522 # elif defined LJ_EWALD_COMB_LB
523 calculate_lj_ewald_comb_LB_F_E(nbparam,
530 # ifdef CALC_ENERGIES
538 # endif /* CALC_ENERGIES */
540 # endif /* LJ_EWALD_COMB_GEOM */
541 # endif /* LJ_EWALD */
543 # ifdef LJ_POT_SWITCH
544 # ifdef CALC_ENERGIES
545 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
547 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
548 # endif /* CALC_ENERGIES */
549 # endif /* LJ_POT_SWITCH */
551 # ifdef VDW_CUTOFF_CHECK
552 /* Separate VDW cut-off check to enable twin-range cut-offs
553 * (rvdw < rcoulomb <= rlist)
555 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
556 F_invr *= vdw_in_range;
557 # ifdef CALC_ENERGIES
558 E_lj_p *= vdw_in_range;
560 # endif /* VDW_CUTOFF_CHECK */
562 # ifdef CALC_ENERGIES
568 # ifdef EXCLUSION_FORCES
569 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
571 F_invr += qi * qj_f * inv_r2 * inv_r;
575 F_invr += qi * qj_f * (int_bit * inv_r2 * inv_r - two_k_rf);
577 # if defined EL_EWALD_ANA
579 * (int_bit * inv_r2 * inv_r + pmecorrF(beta2 * r2) * beta3);
580 # elif defined EL_EWALD_TAB
583 - interpolate_coulomb_force_r(nbparam, r2 * inv_r))
585 # endif /* EL_EWALD_ANA/TAB */
587 # ifdef CALC_ENERGIES
589 E_el += qi * qj_f * (int_bit * inv_r - reactionFieldShift);
593 * (int_bit * inv_r + 0.5f * two_k_rf * r2 - reactionFieldShift);
596 /* 1.0f - erff is faster than erfcf */
598 * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
599 # endif /* EL_EWALD_ANY */
603 /* accumulate j forces in registers */
606 /* accumulate i forces in registers */
611 /* shift the mask bit by 1 */
615 /* reduce j forces */
616 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj, c_fullWarpMask);
620 /* Update the imask with the new one which does not contain the
621 out of range clusters anymore. */
622 pl_cj4[j4].imei[widx].imask = imask;
625 // avoid shared memory WAR hazards between loop iterations
626 __syncwarp(c_fullWarpMask);
629 /* skip central shifts when summing shift forces */
630 if (nb_sci.shift == CENTRAL)
635 float fshift_buf = 0.0f;
637 /* reduce i forces */
638 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
640 ai = (sci * c_nbnxnGpuNumClusterPerSupercluster + i) * c_clSize + tidxi;
641 reduce_force_i_warp_shfl(fci_buf[i], f, &fshift_buf, bCalcFshift, tidxj, ai, c_fullWarpMask);
644 /* add up local shift forces into global mem, tidxj indexes x,y,z */
645 if (bCalcFshift && (tidxj & 3) < 3)
647 float3* fShift = asFloat3(atdat.fShift);
648 atomicAdd(&(fShift[nb_sci.shift].x) + (tidxj & 3), fshift_buf);
651 # ifdef CALC_ENERGIES
652 /* reduce the energies over warps and store into global memory */
653 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx, c_fullWarpMask);
656 #endif /* FUNCTION_DECLARATION_ONLY */
659 #undef MIN_BLOCKS_PER_MP
660 #undef THREADS_PER_BLOCK
663 #undef EXCLUSION_FORCES