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