Introduce ObservablesHistory container
[gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_kernel_fermi.cuh
blobf615978ab2643c981c4940ae73b8da2b09f23d8a
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017, 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  */
36 /*! \internal \file
37  *  \brief
38  *  CUDA non-bonded kernel used through preprocessor-based code generation
39  *  of multiple kernel flavors for CC 2.x, 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  *  \author Berk Hess <hess@kth.se>
45  *  \ingroup module_mdlib
46  */
48 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
49 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/pbcutil/ishift.h"
52 /* Note that floating-point constants in CUDA code should be suffixed
53  * with f (e.g. 0.5f), to stop the compiler producing intermediate
54  * code that is in double precision.
55  */
57 #if GMX_PTX_ARCH >= 300
58 #error "nbnxn_cuda_kernel_fermi.cuh included with GMX_PTX_ARCH >= 300"
59 #endif
61 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
62 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
63 #define EL_EWALD_ANY
64 #endif
66 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
67 /* Macro to control the calculation of exclusion forces in the kernel
68  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
69  * energy terms.
70  *
71  * Note: convenience macro, needs to be undef-ed at the end of the file.
72  */
73 #define EXCLUSION_FORCES
74 #endif
76 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
77 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
78 #define LJ_EWALD
79 #endif
81 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
82 #define LJ_COMB
83 #endif
86    Kernel launch parameters:
87     - #blocks   = #pair lists, blockId = pair list Id
88     - #threads  = c_clSize^2
89     - shmem     = see nbnxn_cuda.cu:calc_shmem_required()
91     Each thread calculates an i force-component taking one pair of i-j atoms.
92  */
94 /**@{*/
95 /*! \brief Definition of kernel launch configuration parameters for CC 2.x.
96  */
98 /* Kernel launch bounds, 16 blocks/multiprocessor can be kept in flight. */
99 #define THREADS_PER_BLOCK   (c_clSize*c_clSize)
101 __launch_bounds__(THREADS_PER_BLOCK)
102 #ifdef PRUNE_NBL
103 #ifdef CALC_ENERGIES
104 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
105 #else
106 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
107 #endif /* CALC_ENERGIES */
108 #else
109 #ifdef CALC_ENERGIES
110 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
111 #else
112 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
113 #endif /* CALC_ENERGIES */
114 #endif /* PRUNE_NBL */
115 (const cu_atomdata_t atdat,
116  const cu_nbparam_t nbparam,
117  const cu_plist_t plist,
118  bool bCalcFshift)
119 #ifdef FUNCTION_DECLARATION_ONLY
120 ;     /* Only do function declaration, omit the function body. */
121 #else
123     /* convenience variables */
124     const nbnxn_sci_t *pl_sci       = plist.sci;
125 #ifndef PRUNE_NBL
126     const
127 #endif
128     nbnxn_cj4_t        *pl_cj4      = plist.cj4;
129     const nbnxn_excl_t *excl        = plist.excl;
130 #ifndef LJ_COMB
131     const int          *atom_types  = atdat.atom_types;
132     int                 ntypes      = atdat.ntypes;
133 #else
134     const float2       *lj_comb     = atdat.lj_comb;
135     float2              ljcp_i, ljcp_j;
136 #endif
137     const float4       *xq          = atdat.xq;
138     float3             *f           = atdat.f;
139     const float3       *shift_vec   = atdat.shift_vec;
140     float               rcoulomb_sq = nbparam.rcoulomb_sq;
141 #ifdef VDW_CUTOFF_CHECK
142     float               rvdw_sq     = nbparam.rvdw_sq;
143     float               vdw_in_range;
144 #endif
145 #ifdef LJ_EWALD
146     float               lje_coeff2, lje_coeff6_6;
147 #endif
148 #ifdef EL_RF
149     float two_k_rf              = nbparam.two_k_rf;
150 #endif
151 #ifdef EL_EWALD_ANA
152     float beta2                 = nbparam.ewald_beta*nbparam.ewald_beta;
153     float beta3                 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
154 #endif
155 #ifdef PRUNE_NBL
156     float rlist_sq              = nbparam.rlist_sq;
157 #endif
159 #ifdef CALC_ENERGIES
160 #ifdef EL_EWALD_ANY
161     float  beta        = nbparam.ewald_beta;
162     float  ewald_shift = nbparam.sh_ewald;
163 #else
164     float  c_rf        = nbparam.c_rf;
165 #endif /* EL_EWALD_ANY */
166     float *e_lj        = atdat.e_lj;
167     float *e_el        = atdat.e_el;
168 #endif /* CALC_ENERGIES */
170     /* thread/block/warp id-s */
171     unsigned int tidxi  = threadIdx.x;
172     unsigned int tidxj  = threadIdx.y;
173     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
174     unsigned int bidx   = blockIdx.x;
175     unsigned int widx   = tidx / warp_size; /* warp index */
177     int          sci, ci, cj,
178                  ai, aj,
179                  cij4_start, cij4_end;
180 #ifndef LJ_COMB
181     int          typei, typej;
182 #endif
183     int          i, jm, j4, wexcl_idx;
184     float        qi, qj_f,
185                  r2, inv_r, inv_r2;
186 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
187     float        inv_r6, c6, c12;
188 #endif
189 #ifdef LJ_COMB_LB
190     float        sigma, epsilon;
191 #endif
192     float        int_bit,
193                  F_invr;
194 #ifdef CALC_ENERGIES
195     float        E_lj, E_el;
196 #endif
197 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
198     float        E_lj_p;
199 #endif
200     unsigned int wexcl, imask, mask_ji;
201     float4       xqbuf;
202     float3       xi, xj, rv, f_ij, fcj_buf;
203     float3       fci_buf[c_numClPerSupercl]; /* i force buffer */
204     nbnxn_sci_t  nb_sci;
206     /*! i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
207     const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
209     /* shmem buffer for i x+q pre-loading */
210     extern __shared__  float4 xqib[];
212     /* shmem buffer for cj, for each warp separately */
213     int   *cjs   = ((int *)(xqib + c_numClPerSupercl * c_clSize));
214     /* shmem j force buffer */
215     float *f_buf = (float *)(cjs + c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize);
217     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
218     sci         = nb_sci.sci;           /* super-cluster */
219     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
220     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
222     {
223         /* Pre-load i-atom x and q into shared memory */
224         ci = sci * c_numClPerSupercl + tidxj;
225         ai = ci * c_clSize + tidxi;
227         xqbuf    = xq[ai] + shift_vec[nb_sci.shift];
228         xqbuf.w *= nbparam.epsfac;
229         xqib[tidxj * c_clSize + tidxi] = xqbuf;
230     }
231     __syncthreads();
233     for (i = 0; i < c_numClPerSupercl; i++)
234     {
235         fci_buf[i] = make_float3(0.0f);
236     }
238 #ifdef LJ_EWALD
239     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
240     lje_coeff2   = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
241     lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*c_oneSixth;
242 #endif
245 #ifdef CALC_ENERGIES
246     E_lj = 0.0f;
247     E_el = 0.0f;
249 #ifdef EXCLUSION_FORCES /* Ewald or RF */
250     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*c_numClPerSupercl)
251     {
252         /* we have the diagonal: add the charge and LJ self interaction energy term */
253         for (i = 0; i < c_numClPerSupercl; i++)
254         {
255 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
256             qi    = xqib[i * c_clSize + tidxi].w;
257             E_el += qi*qi;
258 #endif
260 #ifdef LJ_EWALD
261             E_lj += tex1Dfetch(nbfp_texref, atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2);
262 #endif
263         }
265         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
266 #ifdef LJ_EWALD
267         E_lj /= c_clSize;
268         E_lj *= 0.5f*c_oneSixth*lje_coeff6_6;
269 #endif
271 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
272         /* Correct for epsfac^2 due to adding qi^2 */
273         E_el /= nbparam.epsfac*c_clSize;
274 #if defined EL_RF || defined EL_CUTOFF
275         E_el *= -0.5f*c_rf;
276 #else
277         E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
278 #endif
279 #endif                                  /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
280     }
281 #endif                                  /* EXCLUSION_FORCES */
283 #endif                                  /* CALC_ENERGIES */
285 #ifdef EXCLUSION_FORCES
286     const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
287 #endif
289     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
290     for (j4 = cij4_start; j4 < cij4_end; j4++)
291     {
292         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
293         imask       = pl_cj4[j4].imei[widx].imask;
294         wexcl       = excl[wexcl_idx].pair[(tidx) & (warp_size - 1)];
296 #ifndef PRUNE_NBL
297         if (imask)
298 #endif
299         {
300             /* Pre-load cj into shared memory on both warps separately */
301             if ((tidxj == 0 | tidxj == 4) & (tidxi < c_nbnxnGpuJgroupSize))
302             {
303                 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize/c_splitClSize] = pl_cj4[j4].cj[tidxi];
304             }
306             /* Unrolling this loop with pruning leads to register spilling;
307                Tested with up to nvcc 7.5 */
308 #if !defined PRUNE_NBL
309 #pragma unroll 4
310 #endif
311             for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
312             {
313                 if (imask & (superClInteractionMask << (jm * c_numClPerSupercl)))
314                 {
315                     mask_ji = (1U << (jm * c_numClPerSupercl));
317                     cj      = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize/c_splitClSize];
318                     aj      = cj * c_clSize + tidxj;
320                     /* load j atom data */
321                     xqbuf   = xq[aj];
322                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
323                     qj_f    = xqbuf.w;
324 #ifndef LJ_COMB
325                     typej   = atom_types[aj];
326 #else
327                     ljcp_j  = lj_comb[aj];
328 #endif
330                     fcj_buf = make_float3(0.0f);
332 #if !defined PRUNE_NBL
333 #pragma unroll 8
334 #endif
335                     for (i = 0; i < c_numClPerSupercl; i++)
336                     {
337                         if (imask & mask_ji)
338                         {
339                             ci      = sci * c_numClPerSupercl + i; /* i cluster index */
340                             ai      = ci * c_clSize + tidxi;       /* i atom index */
342                             /* all threads load an atom from i cluster ci into shmem! */
343                             xqbuf   = xqib[i * c_clSize + tidxi];
344                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
346                             /* distance between i and j atoms */
347                             rv      = xi - xj;
348                             r2      = norm2(rv);
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
360                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
362                             /* cutoff & exclusion check */
363 #ifdef EXCLUSION_FORCES
364                             if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
365 #else
366                             if ((r2 < rcoulomb_sq) * int_bit)
367 #endif
368                             {
369                                 /* load the rest of the i-atom parameters */
370                                 qi      = xqbuf.w;
372 #ifndef LJ_COMB
373                                 /* LJ 6*C6 and 12*C12 */
374                                 typei   = atom_types[ai];
375                                 fetch_nbfp_c6_c12(c6, c12, nbparam, ntypes * typei + typej);
376 #else
377                                 ljcp_i  = lj_comb[ai];
378 #ifdef LJ_COMB_GEOM
379                                 c6      = ljcp_i.x * ljcp_j.x;
380                                 c12     = ljcp_i.y * ljcp_j.y;
381 #else
382                                 /* LJ 2^(1/6)*sigma and 12*epsilon */
383                                 sigma   = ljcp_i.x + ljcp_j.x;
384                                 epsilon = ljcp_i.y * ljcp_j.y;
385 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
386                                 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
387 #endif
388 #endif                          /* LJ_COMB_GEOM */
389 #endif                          /* LJ_COMB */
391                                 // Ensure distance do not become so small that r^-12 overflows
392                                 r2      = max(r2, NBNXN_MIN_RSQ);
394                                 inv_r   = rsqrt(r2);
395                                 inv_r2  = inv_r * inv_r;
396 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
397                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
398 #ifdef EXCLUSION_FORCES
399                                 /* We could mask inv_r2, but with Ewald
400                                  * masking both inv_r6 and F_invr is faster */
401                                 inv_r6  *= int_bit;
402 #endif                          /* EXCLUSION_FORCES */
404                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
405 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
406                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*c_oneTwelveth -
407                                                      c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*c_oneSixth);
408 #endif
409 #else                           /* !LJ_COMB_LB || CALC_ENERGIES */
410                                 float sig_r  = sigma*inv_r;
411                                 float sig_r2 = sig_r*sig_r;
412                                 float sig_r6 = sig_r2*sig_r2*sig_r2;
413 #ifdef EXCLUSION_FORCES
414                                 sig_r6 *= int_bit;
415 #endif                          /* EXCLUSION_FORCES */
417                                 F_invr  = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
418 #endif                          /* !LJ_COMB_LB || CALC_ENERGIES */
420 #ifdef LJ_FORCE_SWITCH
421 #ifdef CALC_ENERGIES
422                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
423 #else
424                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
425 #endif /* CALC_ENERGIES */
426 #endif /* LJ_FORCE_SWITCH */
429 #ifdef LJ_EWALD
430 #ifdef LJ_EWALD_COMB_GEOM
431 #ifdef CALC_ENERGIES
432                                 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);
433 #else
434                                 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
435 #endif                          /* CALC_ENERGIES */
436 #elif defined LJ_EWALD_COMB_LB
437                                 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
438 #ifdef CALC_ENERGIES
439                                                                int_bit, &F_invr, &E_lj_p
440 #else
441                                                                0, &F_invr, NULL
442 #endif /* CALC_ENERGIES */
443                                                                );
444 #endif /* LJ_EWALD_COMB_GEOM */
445 #endif /* LJ_EWALD */
447 #ifdef LJ_POT_SWITCH
448 #ifdef CALC_ENERGIES
449                                 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
450 #else
451                                 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
452 #endif /* CALC_ENERGIES */
453 #endif /* LJ_POT_SWITCH */
455 #ifdef VDW_CUTOFF_CHECK
456                                 /* Separate VDW cut-off check to enable twin-range cut-offs
457                                  * (rvdw < rcoulomb <= rlist)
458                                  */
459                                 vdw_in_range  = (r2 < rvdw_sq) ? 1.0f : 0.0f;
460                                 F_invr       *= vdw_in_range;
461 #ifdef CALC_ENERGIES
462                                 E_lj_p       *= vdw_in_range;
463 #endif
464 #endif                          /* VDW_CUTOFF_CHECK */
466 #ifdef CALC_ENERGIES
467                                 E_lj    += E_lj_p;
468 #endif
471 #ifdef EL_CUTOFF
472 #ifdef EXCLUSION_FORCES
473                                 F_invr  += qi * qj_f * int_bit * inv_r2 * inv_r;
474 #else
475                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
476 #endif
477 #endif
478 #ifdef EL_RF
479                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
480 #endif
481 #if defined EL_EWALD_ANA
482                                 F_invr  += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
483 #elif defined EL_EWALD_TAB
484                                 F_invr  += qi * qj_f * (int_bit*inv_r2 -
485                                                         interpolate_coulomb_force_r(nbparam, r2 * inv_r)) * inv_r;
486 #endif                          /* EL_EWALD_ANA/TAB */
488 #ifdef CALC_ENERGIES
489 #ifdef EL_CUTOFF
490                                 E_el    += qi * qj_f * (int_bit*inv_r - c_rf);
491 #endif
492 #ifdef EL_RF
493                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
494 #endif
495 #ifdef EL_EWALD_ANY
496                                 /* 1.0f - erff is faster than erfcf */
497                                 E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
498 #endif                          /* EL_EWALD_ANY */
499 #endif
500                                 f_ij    = rv * F_invr;
502                                 /* accumulate j forces in registers */
503                                 fcj_buf -= f_ij;
505                                 /* accumulate i forces in registers */
506                                 fci_buf[i] += f_ij;
507                             }
508                         }
510                         /* shift the mask bit by 1 */
511                         mask_ji += mask_ji;
512                     }
514                     /* reduce j forces */
515                     /* store j forces in shmem */
516                     f_buf[                   tidx] = fcj_buf.x;
517                     f_buf[    c_fbufStride + tidx] = fcj_buf.y;
518                     f_buf[2 * c_fbufStride + tidx] = fcj_buf.z;
520                     reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
521                 }
522             }
523 #ifdef PRUNE_NBL
524             /* Update the imask with the new one which does not contain the
525                out of range clusters anymore. */
526             pl_cj4[j4].imei[widx].imask = imask;
527 #endif
528         }
529     }
531     /* skip central shifts when summing shift forces */
532     if (nb_sci.shift == CENTRAL)
533     {
534         bCalcFshift = false;
535     }
537     float fshift_buf = 0.0f;
539     /* reduce i forces */
540     for (i = 0; i < c_numClPerSupercl; i++)
541     {
542         ai  = (sci * c_numClPerSupercl + i) * c_clSize + tidxi;
543         f_buf[                   tidx] = fci_buf[i].x;
544         f_buf[    c_fbufStride + tidx] = fci_buf[i].y;
545         f_buf[2 * c_fbufStride + tidx] = fci_buf[i].z;
546         __syncthreads();
547         reduce_force_i(f_buf, f,
548                        &fshift_buf, bCalcFshift,
549                        tidxi, tidxj, ai);
550         __syncthreads();
551     }
553     /* add up local shift forces into global mem, tidxj indexes x,y,z */
554     if (bCalcFshift && tidxj < 3)
555     {
556         atomicAdd(&(atdat.fshift[nb_sci.shift].x) + tidxj, fshift_buf);
557     }
559 #ifdef CALC_ENERGIES
560     /* flush the energies to shmem and reduce them */
561     f_buf[               tidx] = E_lj;
562     f_buf[c_fbufStride + tidx] = E_el;
563     reduce_energy_pow2(f_buf + (tidx & warp_size), e_lj, e_el, tidx & ~warp_size);
564 #endif
566 #endif /* FUNCTION_DECLARATION_ONLY */
568 #undef THREADS_PER_BLOCK
570 #undef EL_EWALD_ANY
571 #undef EXCLUSION_FORCES
572 #undef LJ_EWALD
574 #undef LJ_COMB