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