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, 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/math/utilities.h"
52 #include "gromacs/pbcutil/ishift.h"
53 /* Note that floating-point constants in CUDA code should be suffixed
54 * with f (e.g. 0.5f), to stop the compiler producing intermediate
55 * code that is in double precision.
58 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
59 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
63 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD \
64 || (defined EL_CUTOFF && defined CALC_ENERGIES)
65 /* Macro to control the calculation of exclusion forces in the kernel
66 * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
69 * Note: convenience macro, needs to be undef-ed at the end of the file.
71 # define EXCLUSION_FORCES
74 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
75 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
79 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
84 Kernel launch parameters:
85 - #blocks = #pair lists, blockId = pair list Id
86 - #threads = NTHREAD_Z * c_clSize^2
87 - shmem = see nbnxn_cuda.cu:calc_shmem_required_nonbonded()
89 Each thread calculates an i force-component taking one pair of i-j atoms.
93 /*! \brief Compute capability dependent definition of kernel launch configuration parameters.
95 * NTHREAD_Z controls the number of j-clusters processed concurrently on NTHREAD_Z
96 * warp-pairs per block.
98 * - On CC 3.0-3.5, and >=5.0 NTHREAD_Z == 1, translating to 64 th/block with 16
99 * blocks/multiproc, is the fastest even though this setup gives low occupancy
101 * NTHREAD_Z > 1 results in excessive register spilling unless the minimum blocks
102 * per multiprocessor is reduced proportionally to get the original number of max
103 * threads in flight (and slightly lower performance).
104 * - On CC 3.7 there are enough registers to double the number of threads; using
105 * NTHREADS_Z == 2 is fastest with 16 blocks (TODO: test with RF and other kernels
106 * with low-register use).
108 * Note that the current kernel implementation only supports NTHREAD_Z > 1 with
109 * shuffle-based reduction, hence CC >= 3.0.
112 * NOTEs on Volta / CUDA 9 extensions:
114 * - While active thread masks are required for the warp collectives
115 * (we use any and shfl), the kernel is designed such that all conditions
116 * (other than the inner-most distance check) including loop trip counts
117 * are warp-synchronous. Therefore, we don't need ballot to compute the
118 * active masks as these are all full-warp masks.
120 * - TODO: reconsider the use of __syncwarp(): its only role is currently to prevent
121 * WAR hazard due to the cj preload; we should try to replace it with direct
122 * loads (which may be faster given the improved L1 on Volta).
125 /* Kernel launch bounds for different compute capabilities. The value of NTHREAD_Z
126 * determines the number of threads per block and it is chosen such that
127 * 16 blocks/multiprocessor can be kept in flight.
128 * - CC 3.0,3.5, and >=5.0: NTHREAD_Z=1, (64, 16) bounds
129 * - CC 3.7: NTHREAD_Z=2, (128, 16) bounds
131 * Note: convenience macros, need to be undef-ed at the end of the file.
133 #if GMX_PTX_ARCH == 370
134 # define NTHREAD_Z (2)
135 # define MIN_BLOCKS_PER_MP (16)
137 # define NTHREAD_Z (1)
138 # define MIN_BLOCKS_PER_MP (16)
139 #endif /* GMX_PTX_ARCH == 370 */
140 #define THREADS_PER_BLOCK (c_clSize * c_clSize * NTHREAD_Z)
142 #if GMX_PTX_ARCH >= 350
144 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
146 __launch_bounds__(THREADS_PER_BLOCK)
147 #endif /* GMX_PTX_ARCH >= 350 */
149 # ifdef CALC_ENERGIES
150 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
152 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
153 # endif /* CALC_ENERGIES */
155 # ifdef CALC_ENERGIES
156 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
158 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
159 # endif /* CALC_ENERGIES */
160 #endif /* PRUNE_NBL */
161 (const cu_atomdata_t atdat, const NBParamGpu nbparam, const Nbnxm::gpu_plist plist, bool bCalcFshift)
162 #ifdef FUNCTION_DECLARATION_ONLY
163 ; /* Only do function declaration, omit the function body. */
166 /* convenience variables */
167 const nbnxn_sci_t* pl_sci = plist.sci;
171 nbnxn_cj4_t* pl_cj4 = plist.cj4;
172 const nbnxn_excl_t* excl = plist.excl;
174 const int* atom_types = atdat.atom_types;
175 int ntypes = atdat.ntypes;
177 const float2* lj_comb = atdat.lj_comb;
178 float2 ljcp_i, ljcp_j;
180 const float4* xq = atdat.xq;
182 const float3* shift_vec = atdat.shift_vec;
183 float rcoulomb_sq = nbparam.rcoulomb_sq;
184 # ifdef VDW_CUTOFF_CHECK
185 float rvdw_sq = nbparam.rvdw_sq;
189 float lje_coeff2, lje_coeff6_6;
192 float two_k_rf = nbparam.two_k_rf;
195 float beta2 = nbparam.ewald_beta * nbparam.ewald_beta;
196 float beta3 = nbparam.ewald_beta * nbparam.ewald_beta * nbparam.ewald_beta;
199 float rlist_sq = nbparam.rlistOuter_sq;
202 # ifdef CALC_ENERGIES
204 float beta = nbparam.ewald_beta;
205 float ewald_shift = nbparam.sh_ewald;
207 float c_rf = nbparam.c_rf;
208 # endif /* EL_EWALD_ANY */
209 float* e_lj = atdat.e_lj;
210 float* e_el = atdat.e_el;
211 # endif /* CALC_ENERGIES */
213 /* thread/block/warp id-s */
214 unsigned int tidxi = threadIdx.x;
215 unsigned int tidxj = threadIdx.y;
216 unsigned int tidx = threadIdx.y * blockDim.x + threadIdx.x;
218 unsigned int tidxz = 0;
220 unsigned int tidxz = threadIdx.z;
222 unsigned int bidx = blockIdx.x;
223 unsigned int widx = tidx / warp_size; /* warp index */
225 int sci, ci, cj, ai, aj, cij4_start, cij4_end;
229 int i, jm, j4, wexcl_idx;
230 float qi, qj_f, r2, inv_r, inv_r2;
231 # if !defined LJ_COMB_LB || defined CALC_ENERGIES
232 float inv_r6, c6, c12;
235 float sigma, epsilon;
237 float int_bit, F_invr;
238 # ifdef CALC_ENERGIES
241 # if defined CALC_ENERGIES || defined LJ_POT_SWITCH
244 unsigned int wexcl, imask, mask_ji;
246 float3 xi, xj, rv, f_ij, fcj_buf;
247 float3 fci_buf[c_nbnxnGpuNumClusterPerSupercluster]; /* i force buffer */
250 /*! i-cluster interaction mask for a super-cluster with all c_nbnxnGpuNumClusterPerSupercluster=8 bits set */
251 const unsigned superClInteractionMask = ((1U << c_nbnxnGpuNumClusterPerSupercluster) - 1U);
253 /*********************************************************************
254 * Set up shared memory pointers.
255 * sm_nextSlotPtr should always be updated to point to the "next slot",
256 * that is past the last point where data has been stored.
258 extern __shared__ char sm_dynamicShmem[];
259 char* sm_nextSlotPtr = sm_dynamicShmem;
260 static_assert(sizeof(char) == 1,
261 "The shared memory offset calculation assumes that char is 1 byte");
263 /* shmem buffer for i x+q pre-loading */
264 float4* xqib = (float4*)sm_nextSlotPtr;
265 sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*xqib));
267 /* shmem buffer for cj, for each warp separately */
268 int* cjs = (int*)(sm_nextSlotPtr);
269 /* the cjs buffer's use expects a base pointer offset for pairs of warps in the j-concurrent execution */
270 cjs += tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
271 sm_nextSlotPtr += (NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(*cjs));
274 /* shmem buffer for i atom-type pre-loading */
275 int* atib = (int*)sm_nextSlotPtr;
276 sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*atib));
278 /* shmem buffer for i-atom LJ combination rule parameters */
279 float2* ljcpib = (float2*)sm_nextSlotPtr;
280 sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*ljcpib));
282 /*********************************************************************/
284 nb_sci = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
285 sci = nb_sci.sci; /* super-cluster */
286 cij4_start = nb_sci.cj4_ind_start; /* first ...*/
287 cij4_end = nb_sci.cj4_ind_end; /* and last index of j clusters */
291 /* Pre-load i-atom x and q into shared memory */
292 ci = sci * c_nbnxnGpuNumClusterPerSupercluster + tidxj;
293 ai = ci * c_clSize + tidxi;
295 float* shiftptr = (float*)&shift_vec[nb_sci.shift];
296 xqbuf = xq[ai] + make_float4(LDG(shiftptr), LDG(shiftptr + 1), LDG(shiftptr + 2), 0.0f);
297 xqbuf.w *= nbparam.epsfac;
298 xqib[tidxj * c_clSize + tidxi] = xqbuf;
301 /* Pre-load the i-atom types into shared memory */
302 atib[tidxj * c_clSize + tidxi] = atom_types[ai];
304 /* Pre-load the LJ combination parameters into shared memory */
305 ljcpib[tidxj * c_clSize + tidxi] = lj_comb[ai];
310 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
312 fci_buf[i] = make_float3(0.0f);
316 /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
317 lje_coeff2 = nbparam.ewaldcoeff_lj * nbparam.ewaldcoeff_lj;
318 lje_coeff6_6 = lje_coeff2 * lje_coeff2 * lje_coeff2 * c_oneSixth;
322 # ifdef CALC_ENERGIES
326 # ifdef EXCLUSION_FORCES /* Ewald or RF */
327 if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci * c_nbnxnGpuNumClusterPerSupercluster)
329 /* we have the diagonal: add the charge and LJ self interaction energy term */
330 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
332 # if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
333 qi = xqib[i * c_clSize + tidxi].w;
338 # if DISABLE_CUDA_TEXTURES
339 E_lj += LDG(&nbparam.nbfp[atom_types[(sci * c_nbnxnGpuNumClusterPerSupercluster + i) * c_clSize + tidxi]
340 * (ntypes + 1) * 2]);
342 E_lj += tex1Dfetch<float>(
344 atom_types[(sci * c_nbnxnGpuNumClusterPerSupercluster + i) * c_clSize + tidxi]
350 /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
352 E_lj /= c_clSize * NTHREAD_Z;
353 E_lj *= 0.5f * c_oneSixth * lje_coeff6_6;
356 # if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
357 /* Correct for epsfac^2 due to adding qi^2 */
358 E_el /= nbparam.epsfac * c_clSize * NTHREAD_Z;
359 # if defined EL_RF || defined EL_CUTOFF
360 E_el *= -0.5f * c_rf;
362 E_el *= -beta * M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
364 # endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
366 # endif /* EXCLUSION_FORCES */
368 # endif /* CALC_ENERGIES */
370 # ifdef EXCLUSION_FORCES
371 const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
374 /* loop over the j clusters = seen by any of the atoms in the current super-cluster;
375 * The loop stride NTHREAD_Z ensures that consecutive warps-pairs are assigned
376 * consecutive j4's entries.
378 for (j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
380 wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
381 imask = pl_cj4[j4].imei[widx].imask;
382 wexcl = excl[wexcl_idx].pair[(tidx) & (warp_size - 1)];
388 /* Pre-load cj into shared memory on both warps separately */
389 if ((tidxj == 0 | tidxj == 4) & (tidxi < c_nbnxnGpuJgroupSize))
391 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize / c_splitClSize] = pl_cj4[j4].cj[tidxi];
393 __syncwarp(c_fullWarpMask);
395 /* Unrolling this loop
396 - with pruning leads to register spilling;
397 - on Kepler and later it is much slower;
398 Tested with up to nvcc 7.5 */
399 for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
401 if (imask & (superClInteractionMask << (jm * c_nbnxnGpuNumClusterPerSupercluster)))
403 mask_ji = (1U << (jm * c_nbnxnGpuNumClusterPerSupercluster));
405 cj = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize / c_splitClSize];
406 aj = cj * c_clSize + tidxj;
408 /* load j atom data */
410 xj = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
413 typej = atom_types[aj];
415 ljcp_j = lj_comb[aj];
418 fcj_buf = make_float3(0.0f);
420 # if !defined PRUNE_NBL
423 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
427 ci = sci * c_nbnxnGpuNumClusterPerSupercluster + i; /* i cluster index */
429 /* all threads load an atom from i cluster ci into shmem! */
430 xqbuf = xqib[i * c_clSize + tidxi];
431 xi = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
433 /* distance between i and j atoms */
438 /* If _none_ of the atoms pairs are in cutoff range,
439 the bit corresponding to the current
440 cluster-pair in imask gets set to 0. */
441 if (!__any_sync(c_fullWarpMask, r2 < rlist_sq))
447 int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
449 /* cutoff & exclusion check */
450 # ifdef EXCLUSION_FORCES
451 if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
453 if ((r2 < rcoulomb_sq) * int_bit)
456 /* load the rest of the i-atom parameters */
460 /* LJ 6*C6 and 12*C12 */
461 typei = atib[i * c_clSize + tidxi];
462 fetch_nbfp_c6_c12(c6, c12, nbparam, ntypes * typei + typej);
464 ljcp_i = ljcpib[i * c_clSize + tidxi];
466 c6 = ljcp_i.x * ljcp_j.x;
467 c12 = ljcp_i.y * ljcp_j.y;
469 /* LJ 2^(1/6)*sigma and 12*epsilon */
470 sigma = ljcp_i.x + ljcp_j.x;
471 epsilon = ljcp_i.y * ljcp_j.y;
472 # if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
473 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
475 # endif /* LJ_COMB_GEOM */
476 # endif /* LJ_COMB */
478 // Ensure distance do not become so small that r^-12 overflows
479 r2 = max(r2, c_nbnxnMinDistanceSquared);
482 inv_r2 = inv_r * inv_r;
483 # if !defined LJ_COMB_LB || defined CALC_ENERGIES
484 inv_r6 = inv_r2 * inv_r2 * inv_r2;
485 # ifdef EXCLUSION_FORCES
486 /* We could mask inv_r2, but with Ewald
487 * masking both inv_r6 and F_invr is faster */
489 # endif /* EXCLUSION_FORCES */
491 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
492 # if defined CALC_ENERGIES || defined LJ_POT_SWITCH
494 * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot) * c_oneTwelveth
495 - c6 * (inv_r6 + nbparam.dispersion_shift.cpot) * c_oneSixth);
497 # else /* !LJ_COMB_LB || CALC_ENERGIES */
498 float sig_r = sigma * inv_r;
499 float sig_r2 = sig_r * sig_r;
500 float sig_r6 = sig_r2 * sig_r2 * sig_r2;
501 # ifdef EXCLUSION_FORCES
503 # endif /* EXCLUSION_FORCES */
505 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
506 # endif /* !LJ_COMB_LB || CALC_ENERGIES */
508 # ifdef LJ_FORCE_SWITCH
509 # ifdef CALC_ENERGIES
510 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
512 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
513 # endif /* CALC_ENERGIES */
514 # endif /* LJ_FORCE_SWITCH */
518 # ifdef LJ_EWALD_COMB_GEOM
519 # ifdef CALC_ENERGIES
520 calculate_lj_ewald_comb_geom_F_E(nbparam, typei, typej, r2, inv_r2,
521 lje_coeff2, lje_coeff6_6, int_bit,
524 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2,
525 lje_coeff2, lje_coeff6_6, &F_invr);
526 # endif /* CALC_ENERGIES */
527 # elif defined LJ_EWALD_COMB_LB
528 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2,
529 lje_coeff2, lje_coeff6_6,
530 # ifdef CALC_ENERGIES
531 int_bit, &F_invr, &E_lj_p
534 # endif /* CALC_ENERGIES */
536 # endif /* LJ_EWALD_COMB_GEOM */
537 # endif /* LJ_EWALD */
539 # ifdef LJ_POT_SWITCH
540 # ifdef CALC_ENERGIES
541 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
543 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
544 # endif /* CALC_ENERGIES */
545 # endif /* LJ_POT_SWITCH */
547 # ifdef VDW_CUTOFF_CHECK
548 /* Separate VDW cut-off check to enable twin-range cut-offs
549 * (rvdw < rcoulomb <= rlist)
551 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
552 F_invr *= vdw_in_range;
553 # ifdef CALC_ENERGIES
554 E_lj_p *= vdw_in_range;
556 # endif /* VDW_CUTOFF_CHECK */
558 # ifdef CALC_ENERGIES
564 # ifdef EXCLUSION_FORCES
565 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
567 F_invr += qi * qj_f * inv_r2 * inv_r;
571 F_invr += qi * qj_f * (int_bit * inv_r2 * inv_r - two_k_rf);
573 # if defined EL_EWALD_ANA
575 * (int_bit * inv_r2 * inv_r + pmecorrF(beta2 * r2) * beta3);
576 # elif defined EL_EWALD_TAB
579 - interpolate_coulomb_force_r(nbparam, r2 * inv_r))
581 # endif /* EL_EWALD_ANA/TAB */
583 # ifdef CALC_ENERGIES
585 E_el += qi * qj_f * (int_bit * inv_r - c_rf);
588 E_el += qi * qj_f * (int_bit * inv_r + 0.5f * two_k_rf * r2 - c_rf);
591 /* 1.0f - erff is faster than erfcf */
593 * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
594 # endif /* EL_EWALD_ANY */
598 /* accumulate j forces in registers */
601 /* accumulate i forces in registers */
606 /* shift the mask bit by 1 */
610 /* reduce j forces */
611 reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj, c_fullWarpMask);
615 /* Update the imask with the new one which does not contain the
616 out of range clusters anymore. */
617 pl_cj4[j4].imei[widx].imask = imask;
620 // avoid shared memory WAR hazards between loop iterations
621 __syncwarp(c_fullWarpMask);
624 /* skip central shifts when summing shift forces */
625 if (nb_sci.shift == CENTRAL)
630 float fshift_buf = 0.0f;
632 /* reduce i forces */
633 for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
635 ai = (sci * c_nbnxnGpuNumClusterPerSupercluster + i) * c_clSize + tidxi;
636 reduce_force_i_warp_shfl(fci_buf[i], f, &fshift_buf, bCalcFshift, tidxj, ai, c_fullWarpMask);
639 /* add up local shift forces into global mem, tidxj indexes x,y,z */
640 if (bCalcFshift && (tidxj & 3) < 3)
642 atomicAdd(&(atdat.fshift[nb_sci.shift].x) + (tidxj & 3), fshift_buf);
645 # ifdef CALC_ENERGIES
646 /* reduce the energies over warps and store into global memory */
647 reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx, c_fullWarpMask);
650 #endif /* FUNCTION_DECLARATION_ONLY */
653 #undef MIN_BLOCKS_PER_MP
654 #undef THREADS_PER_BLOCK
657 #undef EXCLUSION_FORCES