Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_kernel.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 #include "gromacs/math/utilities.h"
37 #include "gromacs/pbcutil/ishift.h"
38 /* Note that floating-point constants in CUDA code should be suffixed
39  * with f (e.g. 0.5f), to stop the compiler producing intermediate
40  * code that is in double precision.
41  */
42
43 #if __CUDA_ARCH__ >= 300
44 #define REDUCE_SHUFFLE
45 /* On Kepler pre-loading i-atom types to shmem gives a few %,
46    but on Fermi it does not */
47 #define IATYPE_SHMEM
48 #endif
49
50 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
51 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
52 #define EL_EWALD_ANY
53 #endif
54
55 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
56 /* Macro to control the calculation of exclusion forces in the kernel
57  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
58  * energy terms.
59  *
60  * Note: convenience macro, needs to be undef-ed at the end of the file.
61  */
62 #define EXCLUSION_FORCES
63 #endif
64
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. */
67 #define LJ_EWALD
68 #endif
69
70 /*
71    Kernel launch parameters:
72     - #blocks   = #pair lists, blockId = pair list Id
73     - #threads  = CL_SIZE^2
74     - shmem     = CL_SIZE^2 * sizeof(float)
75
76     Each thread calculates an i force-component taking one pair of i-j atoms.
77  */
78 #if __CUDA_ARCH__ >= 350
79 __launch_bounds__(64, 16)
80 #endif
81 #ifdef PRUNE_NBL
82 #ifdef CALC_ENERGIES
83 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
84 #else
85 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
86 #endif
87 #else
88 #ifdef CALC_ENERGIES
89 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
90 #else
91 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
92 #endif
93 #endif
94 (const cu_atomdata_t atdat,
95  const cu_nbparam_t nbparam,
96  const cu_plist_t plist,
97  bool bCalcFshift)
98 {
99     /* convenience variables */
100     const nbnxn_sci_t *pl_sci       = plist.sci;
101 #ifndef PRUNE_NBL
102     const
103 #endif
104     nbnxn_cj4_t        *pl_cj4      = plist.cj4;
105     const nbnxn_excl_t *excl        = plist.excl;
106     const int          *atom_types  = atdat.atom_types;
107     int                 ntypes      = atdat.ntypes;
108     const float4       *xq          = atdat.xq;
109     float3             *f           = atdat.f;
110     const float3       *shift_vec   = atdat.shift_vec;
111     float               rcoulomb_sq = nbparam.rcoulomb_sq;
112 #ifdef VDW_CUTOFF_CHECK
113     float               rvdw_sq     = nbparam.rvdw_sq;
114     float               vdw_in_range;
115 #endif
116 #ifdef LJ_EWALD
117     float               lje_coeff2, lje_coeff6_6;
118 #endif
119 #ifdef EL_RF
120     float two_k_rf              = nbparam.two_k_rf;
121 #endif
122 #ifdef EL_EWALD_TAB
123     float coulomb_tab_scale     = nbparam.coulomb_tab_scale;
124 #endif
125 #ifdef EL_EWALD_ANA
126     float beta2                 = nbparam.ewald_beta*nbparam.ewald_beta;
127     float beta3                 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
128 #endif
129 #ifdef PRUNE_NBL
130     float rlist_sq              = nbparam.rlist_sq;
131 #endif
132
133 #ifdef CALC_ENERGIES
134 #ifdef EL_EWALD_ANY
135     float  beta        = nbparam.ewald_beta;
136     float  ewald_shift = nbparam.sh_ewald;
137 #else
138     float  c_rf        = nbparam.c_rf;
139 #endif /* EL_EWALD_ANY */
140     float *e_lj        = atdat.e_lj;
141     float *e_el        = atdat.e_el;
142 #endif /* CALC_ENERGIES */
143
144     /* thread/block/warp id-s */
145     unsigned int tidxi  = threadIdx.x;
146     unsigned int tidxj  = threadIdx.y;
147     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
148     unsigned int bidx   = blockIdx.x;
149     unsigned int widx   = tidx / WARP_SIZE; /* warp index */
150
151     int          sci, ci, cj, ci_offset,
152                  ai, aj,
153                  cij4_start, cij4_end,
154                  typei, typej,
155                  i, jm, j4, wexcl_idx;
156     float        qi, qj_f,
157                  r2, inv_r, inv_r2, inv_r6,
158                  c6, c12,
159                  int_bit,
160                  F_invr;
161 #ifdef CALC_ENERGIES
162     float        E_lj, E_el;
163 #endif
164 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
165     float        E_lj_p;
166 #endif
167     unsigned int wexcl, imask, mask_ji;
168     float4       xqbuf;
169     float3       xi, xj, rv, f_ij, fcj_buf, fshift_buf;
170     float3       fci_buf[NCL_PER_SUPERCL]; /* i force buffer */
171     nbnxn_sci_t  nb_sci;
172
173     /* shmem buffer for i x+q pre-loading */
174     extern __shared__  float4 xqib[];
175     /* shmem buffer for cj, for both warps separately */
176     int *cjs     = (int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
177 #ifdef IATYPE_SHMEM
178     /* shmem buffer for i atom-type pre-loading */
179     int *atib = (int *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
180 #endif
181
182 #ifndef REDUCE_SHUFFLE
183     /* shmem j force buffer */
184 #ifdef IATYPE_SHMEM
185     float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
186 #else
187     float *f_buf = (float *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
188 #endif
189 #endif
190
191     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
192     sci         = nb_sci.sci;           /* super-cluster */
193     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
194     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
195
196     /* Pre-load i-atom x and q into shared memory */
197     ci = sci * NCL_PER_SUPERCL + tidxj;
198     ai = ci * CL_SIZE + tidxi;
199     xqib[tidxj * CL_SIZE + tidxi] = xq[ai] + shift_vec[nb_sci.shift];
200 #ifdef IATYPE_SHMEM
201     /* Pre-load the i-atom types into shared memory */
202     atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
203 #endif
204     __syncthreads();
205
206     for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
207     {
208         fci_buf[ci_offset] = make_float3(0.0f);
209     }
210
211 #ifdef LJ_EWALD
212     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
213     lje_coeff2   = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
214     lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
215 #endif /* LJ_EWALD */
216
217
218 #ifdef CALC_ENERGIES
219     E_lj = 0.0f;
220     E_el = 0.0f;
221
222 #if defined EXCLUSION_FORCES /* Ewald or RF */
223     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
224     {
225         /* we have the diagonal: add the charge and LJ self interaction energy term */
226         for (i = 0; i < NCL_PER_SUPERCL; i++)
227         {
228 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
229             qi    = xqib[i * CL_SIZE + tidxi].w;
230             E_el += qi*qi;
231 #endif
232
233 #if defined LJ_EWALD
234 #ifdef USE_TEXOBJ
235             E_lj += tex1Dfetch<float>(nbparam.nbfp_texobj, atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2);
236 #else
237             E_lj += tex1Dfetch(nbfp_texref, atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2);
238 #endif /* USE_TEXOBJ */
239 #endif /* LJ_EWALD */
240
241         }
242
243         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
244 #ifdef LJ_EWALD
245         E_lj /= CL_SIZE;
246         E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
247 #endif  /* LJ_EWALD */
248
249 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
250         E_el /= CL_SIZE;
251 #if defined EL_RF || defined EL_CUTOFF
252         E_el *= -nbparam.epsfac*0.5f*c_rf;
253 #else
254         E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
255 #endif
256 #endif                                                 /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
257     }
258 #endif                                                 /* EXCLUSION_FORCES */
259
260 #endif                                                 /* CALC_ENERGIES */
261
262     /* skip central shifts when summing shift forces */
263     if (nb_sci.shift == CENTRAL)
264     {
265         bCalcFshift = false;
266     }
267
268     fshift_buf = make_float3(0.0f);
269
270     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
271     for (j4 = cij4_start; j4 < cij4_end; j4++)
272     {
273         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
274         imask       = pl_cj4[j4].imei[widx].imask;
275         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
276
277 #ifndef PRUNE_NBL
278         if (imask)
279 #endif
280         {
281             /* Pre-load cj into shared memory on both warps separately */
282             if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
283             {
284                 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
285             }
286
287             /* Unrolling this loop
288                - with pruning leads to register spilling;
289                - on Kepler is much slower;
290                - doesn't work on CUDA <v4.1
291                Tested with nvcc 3.2 - 5.0.7 */
292 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && CUDA_VERSION >= 4010
293 #pragma unroll 4
294 #endif
295             for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
296             {
297                 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
298                 {
299                     mask_ji = (1U << (jm * NCL_PER_SUPERCL));
300
301                     cj      = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
302                     aj      = cj * CL_SIZE + tidxj;
303
304                     /* load j atom data */
305                     xqbuf   = xq[aj];
306                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
307                     qj_f    = nbparam.epsfac * xqbuf.w;
308                     typej   = atom_types[aj];
309
310                     fcj_buf = make_float3(0.0f);
311
312                     /* The PME and RF kernels don't unroll with CUDA <v4.1. */
313 #if !defined PRUNE_NBL && !(CUDA_VERSION < 4010 && defined EXCLUSION_FORCES)
314 #pragma unroll 8
315 #endif
316                     for (i = 0; i < NCL_PER_SUPERCL; i++)
317                     {
318                         if (imask & mask_ji)
319                         {
320                             ci_offset   = i;                     /* i force buffer offset */
321
322                             ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
323                             ai      = ci * CL_SIZE + tidxi;      /* i atom index */
324
325                             /* all threads load an atom from i cluster ci into shmem! */
326                             xqbuf   = xqib[i * CL_SIZE + tidxi];
327                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
328
329                             /* distance between i and j atoms */
330                             rv      = xi - xj;
331                             r2      = norm2(rv);
332
333 #ifdef PRUNE_NBL
334                             /* If _none_ of the atoms pairs are in cutoff range,
335                                the bit corresponding to the current
336                                cluster-pair in imask gets set to 0. */
337                             if (!__any(r2 < rlist_sq))
338                             {
339                                 imask &= ~mask_ji;
340                             }
341 #endif
342
343                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
344
345                             /* cutoff & exclusion check */
346 #ifdef EXCLUSION_FORCES
347                             if (r2 < rcoulomb_sq *
348                                 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
349 #else
350                             if (r2 < rcoulomb_sq * int_bit)
351 #endif
352                             {
353                                 /* load the rest of the i-atom parameters */
354                                 qi      = xqbuf.w;
355 #ifdef IATYPE_SHMEM
356                                 typei   = atib[i * CL_SIZE + tidxi];
357 #else
358                                 typei   = atom_types[ai];
359 #endif
360
361                                 /* LJ 6*C6 and 12*C12 */
362 #ifdef USE_TEXOBJ
363                                 c6      = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej));
364                                 c12     = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej) + 1);
365 #else
366                                 c6      = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
367                                 c12     = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
368 #endif                          /* USE_TEXOBJ */
369
370
371                                 /* avoid NaN for excluded pairs at r=0 */
372                                 r2      += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
373
374                                 inv_r   = rsqrt(r2);
375                                 inv_r2  = inv_r * inv_r;
376                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
377 #if defined EXCLUSION_FORCES
378                                 /* We could mask inv_r2, but with Ewald
379                                  * masking both inv_r6 and F_invr is faster */
380                                 inv_r6  *= int_bit;
381 #endif                          /* EXCLUSION_FORCES */
382
383                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
384 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
385                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*ONE_TWELVETH_F -
386                                                      c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*ONE_SIXTH_F);
387 #endif
388
389 #ifdef LJ_FORCE_SWITCH
390 #ifdef CALC_ENERGIES
391                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
392 #else
393                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
394 #endif /* CALC_ENERGIES */
395 #endif /* LJ_FORCE_SWITCH */
396
397
398 #ifdef LJ_EWALD
399 #ifdef LJ_EWALD_COMB_GEOM
400 #ifdef CALC_ENERGIES
401                                 calculate_lj_ewald_comb_geom_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
402 #else
403                                 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
404 #endif                          /* CALC_ENERGIES */
405 #elif defined LJ_EWALD_COMB_LB
406                                 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
407 #ifdef CALC_ENERGIES
408                                                                int_bit, &F_invr, &E_lj_p
409 #else
410                                                                0, &F_invr, NULL
411 #endif /* CALC_ENERGIES */
412                                                                );
413 #endif /* LJ_EWALD_COMB_GEOM */
414 #endif /* LJ_EWALD */
415
416 #ifdef VDW_CUTOFF_CHECK
417                                 /* Separate VDW cut-off check to enable twin-range cut-offs
418                                  * (rvdw < rcoulomb <= rlist)
419                                  */
420                                 vdw_in_range  = (r2 < rvdw_sq) ? 1.0f : 0.0f;
421                                 F_invr       *= vdw_in_range;
422 #ifdef CALC_ENERGIES
423                                 E_lj_p       *= vdw_in_range;
424 #endif
425 #endif                          /* VDW_CUTOFF_CHECK */
426
427 #ifdef LJ_POT_SWITCH
428 #ifdef CALC_ENERGIES
429                                 calculate_potential_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
430 #else
431                                 calculate_potential_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
432 #endif /* CALC_ENERGIES */
433 #endif /* LJ_POT_SWITCH */
434
435 #ifdef CALC_ENERGIES
436                                 E_lj    += E_lj_p;
437 #endif
438
439
440 #ifdef EL_CUTOFF
441 #ifdef EXCLUSION_FORCES
442                                 F_invr  += qi * qj_f * int_bit * inv_r2 * inv_r;
443 #else
444                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
445 #endif
446 #endif
447 #ifdef EL_RF
448                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
449 #endif
450 #if defined EL_EWALD_ANA
451                                 F_invr  += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
452 #elif defined EL_EWALD_TAB
453                                 F_invr  += qi * qj_f * (int_bit*inv_r2 -
454 #ifdef USE_TEXOBJ
455                                                         interpolate_coulomb_force_r(nbparam.coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
456 #else
457                                                         interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)
458 #endif /* USE_TEXOBJ */
459                                                         ) * inv_r;
460 #endif /* EL_EWALD_ANA/TAB */
461
462 #ifdef CALC_ENERGIES
463 #ifdef EL_CUTOFF
464                                 E_el    += qi * qj_f * (int_bit*inv_r - c_rf);
465 #endif
466 #ifdef EL_RF
467                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
468 #endif
469 #ifdef EL_EWALD_ANY
470                                 /* 1.0f - erff is faster than erfcf */
471                                 E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
472 #endif                          /* EL_EWALD_ANY */
473 #endif
474                                 f_ij    = rv * F_invr;
475
476                                 /* accumulate j forces in registers */
477                                 fcj_buf -= f_ij;
478
479                                 /* accumulate i forces in registers */
480                                 fci_buf[ci_offset] += f_ij;
481                             }
482                         }
483
484                         /* shift the mask bit by 1 */
485                         mask_ji += mask_ji;
486                     }
487
488                     /* reduce j forces */
489 #ifdef REDUCE_SHUFFLE
490                     reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
491 #else
492                     /* store j forces in shmem */
493                     f_buf[                  tidx] = fcj_buf.x;
494                     f_buf[    FBUF_STRIDE + tidx] = fcj_buf.y;
495                     f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
496
497                     reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
498 #endif
499                 }
500             }
501 #ifdef PRUNE_NBL
502             /* Update the imask with the new one which does not contain the
503                out of range clusters anymore. */
504             pl_cj4[j4].imei[widx].imask = imask;
505 #endif
506         }
507     }
508
509     /* reduce i forces */
510     for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
511     {
512         ai  = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
513 #ifdef REDUCE_SHUFFLE
514         reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
515                                  &fshift_buf, bCalcFshift,
516                                  tidxj, ai);
517 #else
518         f_buf[                  tidx] = fci_buf[ci_offset].x;
519         f_buf[    FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
520         f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
521         __syncthreads();
522         reduce_force_i(f_buf, f,
523                        &fshift_buf, bCalcFshift,
524                        tidxi, tidxj, ai);
525         __syncthreads();
526 #endif
527     }
528
529     /* add up local shift forces into global mem */
530 #ifdef REDUCE_SHUFFLE
531     if (bCalcFshift && (tidxj == 0 || tidxj == 4))
532 #else
533     if (bCalcFshift && tidxj == 0)
534 #endif
535     {
536         atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
537         atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
538         atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
539     }
540
541 #ifdef CALC_ENERGIES
542 #ifdef REDUCE_SHUFFLE
543     /* reduce the energies over warps and store into global memory */
544     reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
545 #else
546     /* flush the energies to shmem and reduce them */
547     f_buf[              tidx] = E_lj;
548     f_buf[FBUF_STRIDE + tidx] = E_el;
549     reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
550 #endif
551 #endif
552 }
553
554 #undef EL_EWALD_ANY
555 #undef EXCLUSION_FORCES
556 #undef LJ_EWALD