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