Remove variable precision gro writing
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_kernel_fermi.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016, 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 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  */
47
48 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/pbcutil/ishift.h"
51 /* Note that floating-point constants in CUDA code should be suffixed
52  * with f (e.g. 0.5f), to stop the compiler producing intermediate
53  * code that is in double precision.
54  */
55
56 #if GMX_PTX_ARCH >= 300
57 #error "nbnxn_cuda_kernel_fermi.cuh included with GMX_PTX_ARCH >= 300"
58 #endif
59
60 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
61 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
62 #define EL_EWALD_ANY
63 #endif
64
65 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
66 /* Macro to control the calculation of exclusion forces in the kernel
67  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
68  * energy terms.
69  *
70  * Note: convenience macro, needs to be undef-ed at the end of the file.
71  */
72 #define EXCLUSION_FORCES
73 #endif
74
75 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
76 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
77 #define LJ_EWALD
78 #endif
79
80 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
81 #define LJ_COMB
82 #endif
83
84 /*
85    Kernel launch parameters:
86     - #blocks   = #pair lists, blockId = pair list Id
87     - #threads  = c_clSize^2
88     - shmem     = see nbnxn_cuda.cu:calc_shmem_required()
89
90     Each thread calculates an i force-component taking one pair of i-j atoms.
91  */
92
93 /**@{*/
94 /*! \brief Definition of kernel launch configuration parameters for CC 2.x.
95  */
96
97 /* Kernel launch bounds, 16 blocks/multiprocessor can be kept in flight. */
98 #define THREADS_PER_BLOCK   (c_clSize*c_clSize)
99
100 __launch_bounds__(THREADS_PER_BLOCK)
101 #ifdef PRUNE_NBL
102 #ifdef CALC_ENERGIES
103 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
104 #else
105 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
106 #endif /* CALC_ENERGIES */
107 #else
108 #ifdef CALC_ENERGIES
109 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
110 #else
111 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
112 #endif /* CALC_ENERGIES */
113 #endif /* PRUNE_NBL */
114 (const cu_atomdata_t atdat,
115  const cu_nbparam_t nbparam,
116  const cu_plist_t plist,
117  bool bCalcFshift)
118 #ifdef FUNCTION_DECLARATION_ONLY
119 ;     /* Only do function declaration, omit the function body. */
120 #else
121 {
122     /* convenience variables */
123     const nbnxn_sci_t *pl_sci       = plist.sci;
124 #ifndef PRUNE_NBL
125     const
126 #endif
127     nbnxn_cj4_t        *pl_cj4      = plist.cj4;
128     const nbnxn_excl_t *excl        = plist.excl;
129 #ifndef LJ_COMB
130     const int          *atom_types  = atdat.atom_types;
131     int                 ntypes      = atdat.ntypes;
132 #else
133     const float2       *lj_comb     = atdat.lj_comb;
134     float2              ljcp_i, ljcp_j;
135 #endif
136     const float4       *xq          = atdat.xq;
137     float3             *f           = atdat.f;
138     const float3       *shift_vec   = atdat.shift_vec;
139     float               rcoulomb_sq = nbparam.rcoulomb_sq;
140 #ifdef VDW_CUTOFF_CHECK
141     float               rvdw_sq     = nbparam.rvdw_sq;
142     float               vdw_in_range;
143 #endif
144 #ifdef LJ_EWALD
145     float               lje_coeff2, lje_coeff6_6;
146 #endif
147 #ifdef EL_RF
148     float two_k_rf              = nbparam.two_k_rf;
149 #endif
150 #ifdef EL_EWALD_TAB
151     float coulomb_tab_scale     = nbparam.coulomb_tab_scale;
152 #endif
153 #ifdef EL_EWALD_ANA
154     float beta2                 = nbparam.ewald_beta*nbparam.ewald_beta;
155     float beta3                 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
156 #endif
157 #ifdef PRUNE_NBL
158     float rlist_sq              = nbparam.rlist_sq;
159 #endif
160
161 #ifdef CALC_ENERGIES
162 #ifdef EL_EWALD_ANY
163     float  beta        = nbparam.ewald_beta;
164     float  ewald_shift = nbparam.sh_ewald;
165 #else
166     float  c_rf        = nbparam.c_rf;
167 #endif /* EL_EWALD_ANY */
168     float *e_lj        = atdat.e_lj;
169     float *e_el        = atdat.e_el;
170 #endif /* CALC_ENERGIES */
171
172     /* thread/block/warp id-s */
173     unsigned int tidxi  = threadIdx.x;
174     unsigned int tidxj  = threadIdx.y;
175     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
176     unsigned int bidx   = blockIdx.x;
177     unsigned int widx   = tidx / warp_size; /* warp index */
178
179     int          sci, ci, cj,
180                  ai, aj,
181                  cij4_start, cij4_end;
182 #ifndef LJ_COMB
183     int          typei, typej;
184 #endif
185     int          i, jm, j4, wexcl_idx;
186     float        qi, qj_f,
187                  r2, inv_r, inv_r2;
188 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
189     float        inv_r6, c6, c12;
190 #endif
191 #ifdef LJ_COMB_LB
192     float        sigma, epsilon;
193 #endif
194     float        int_bit,
195                  F_invr;
196 #ifdef CALC_ENERGIES
197     float        E_lj, E_el;
198 #endif
199 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
200     float        E_lj_p;
201 #endif
202     unsigned int wexcl, imask, mask_ji;
203     float4       xqbuf;
204     float3       xi, xj, rv, f_ij, fcj_buf;
205     float3       fci_buf[c_numClPerSupercl]; /* i force buffer */
206     nbnxn_sci_t  nb_sci;
207
208     /*! i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
209     const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
210
211     /* shmem buffer for i x+q pre-loading */
212     extern __shared__  float4 xqib[];
213
214     /* shmem buffer for cj, for each warp separately */
215     int   *cjs   = ((int *)(xqib + c_numClPerSupercl * c_clSize));
216     /* shmem j force buffer */
217     float *f_buf = (float *)(cjs + c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize);
218
219     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
220     sci         = nb_sci.sci;           /* super-cluster */
221     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
222     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
223
224     {
225         /* Pre-load i-atom x and q into shared memory */
226         ci = sci * c_numClPerSupercl + tidxj;
227         ai = ci * c_clSize + tidxi;
228
229         xqbuf    = xq[ai] + shift_vec[nb_sci.shift];
230         xqbuf.w *= nbparam.epsfac;
231         xqib[tidxj * c_clSize + tidxi] = xqbuf;
232     }
233     __syncthreads();
234
235     for (i = 0; i < c_numClPerSupercl; i++)
236     {
237         fci_buf[i] = make_float3(0.0f);
238     }
239
240 #ifdef LJ_EWALD
241     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
242     lje_coeff2   = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
243     lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*c_oneSixth;
244 #endif
245
246
247 #ifdef CALC_ENERGIES
248     E_lj = 0.0f;
249     E_el = 0.0f;
250
251 #ifdef EXCLUSION_FORCES /* Ewald or RF */
252     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*c_numClPerSupercl)
253     {
254         /* we have the diagonal: add the charge and LJ self interaction energy term */
255         for (i = 0; i < c_numClPerSupercl; i++)
256         {
257 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
258             qi    = xqib[i * c_clSize + tidxi].w;
259             E_el += qi*qi;
260 #endif
261
262 #ifdef LJ_EWALD
263             E_lj += tex1Dfetch(nbfp_texref, atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2);
264 #endif
265         }
266
267         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
268 #ifdef LJ_EWALD
269         E_lj /= c_clSize;
270         E_lj *= 0.5f*c_oneSixth*lje_coeff6_6;
271 #endif
272
273 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
274         /* Correct for epsfac^2 due to adding qi^2 */
275         E_el /= nbparam.epsfac*c_clSize;
276 #if defined EL_RF || defined EL_CUTOFF
277         E_el *= -0.5f*c_rf;
278 #else
279         E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
280 #endif
281 #endif                                  /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
282     }
283 #endif                                  /* EXCLUSION_FORCES */
284
285 #endif                                  /* CALC_ENERGIES */
286
287     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
288     for (j4 = cij4_start; j4 < cij4_end; j4++)
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 < c_nbnxnGpuJgroupSize)
300             {
301                 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize/c_splitClSize] = pl_cj4[j4].cj[tidxi];
302             }
303
304             /* Unrolling this loop with pruning leads to register spilling;
305                Tested with up to nvcc 7.5 */
306 #if !defined PRUNE_NBL
307 #pragma unroll 4
308 #endif
309             for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
310             {
311                 if (imask & (superClInteractionMask << (jm * c_numClPerSupercl)))
312                 {
313                     mask_ji = (1U << (jm * c_numClPerSupercl));
314
315                     cj      = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize/c_splitClSize];
316                     aj      = cj * c_clSize + tidxj;
317
318                     /* load j atom data */
319                     xqbuf   = xq[aj];
320                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
321                     qj_f    = xqbuf.w;
322 #ifndef LJ_COMB
323                     typej   = atom_types[aj];
324 #else
325                     ljcp_j  = lj_comb[aj];
326 #endif
327
328                     fcj_buf = make_float3(0.0f);
329
330 #if !defined PRUNE_NBL
331 #pragma unroll 8
332 #endif
333                     for (i = 0; i < c_numClPerSupercl; i++)
334                     {
335                         if (imask & mask_ji)
336                         {
337                             ci      = sci * c_numClPerSupercl + i; /* i cluster index */
338                             ai      = ci * c_clSize + tidxi;       /* i atom index */
339
340                             /* all threads load an atom from i cluster ci into shmem! */
341                             xqbuf   = xqib[i * c_clSize + tidxi];
342                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
343
344                             /* distance between i and j atoms */
345                             rv      = xi - xj;
346                             r2      = norm2(rv);
347
348 #ifdef PRUNE_NBL
349                             /* If _none_ of the atoms pairs are in cutoff range,
350                                the bit corresponding to the current
351                                cluster-pair in imask gets set to 0. */
352                             if (!__any(r2 < rlist_sq))
353                             {
354                                 imask &= ~mask_ji;
355                             }
356 #endif
357
358                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
359
360                             /* cutoff & exclusion check */
361 #ifdef EXCLUSION_FORCES
362                             if (r2 < rcoulomb_sq *
363                                 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
364 #else
365                             if (r2 < rcoulomb_sq * int_bit)
366 #endif
367                             {
368                                 /* load the rest of the i-atom parameters */
369                                 qi      = xqbuf.w;
370
371 #ifndef LJ_COMB
372                                 /* LJ 6*C6 and 12*C12 */
373                                 typei   = atom_types[ai];
374                                 c6      = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
375                                 c12     = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
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 */
390
391                                 // Ensure distance do not become so small that r^-12 overflows
392                                 r2      = max(r2, NBNXN_MIN_RSQ);
393
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 */
403
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 */
416
417                                 F_invr  = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
418 #endif                          /* !LJ_COMB_LB || CALC_ENERGIES */
419
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 */
427
428
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 */
446
447 #ifdef VDW_CUTOFF_CHECK
448                                 /* Separate VDW cut-off check to enable twin-range cut-offs
449                                  * (rvdw < rcoulomb <= rlist)
450                                  */
451                                 vdw_in_range  = (r2 < rvdw_sq) ? 1.0f : 0.0f;
452                                 F_invr       *= vdw_in_range;
453 #ifdef CALC_ENERGIES
454                                 E_lj_p       *= vdw_in_range;
455 #endif
456 #endif                          /* VDW_CUTOFF_CHECK */
457
458 #ifdef LJ_POT_SWITCH
459 #ifdef CALC_ENERGIES
460                                 calculate_potential_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
461 #else
462                                 calculate_potential_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
463 #endif /* CALC_ENERGIES */
464 #endif /* LJ_POT_SWITCH */
465
466 #ifdef CALC_ENERGIES
467                                 E_lj    += E_lj_p;
468 #endif
469
470
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(r2 * inv_r, coulomb_tab_scale)) * inv_r;
486 #endif                          /* EL_EWALD_ANA/TAB */
487
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;
501
502                                 /* accumulate j forces in registers */
503                                 fcj_buf -= f_ij;
504
505                                 /* accumulate i forces in registers */
506                                 fci_buf[i] += f_ij;
507                             }
508                         }
509
510                         /* shift the mask bit by 1 */
511                         mask_ji += mask_ji;
512                     }
513
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;
519
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     }
530
531     /* skip central shifts when summing shift forces */
532     if (nb_sci.shift == CENTRAL)
533     {
534         bCalcFshift = false;
535     }
536
537     float fshift_buf = 0.0f;
538
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     }
552
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     }
558
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
565 }
566 #endif /* FUNCTION_DECLARATION_ONLY */
567
568 #undef THREADS_PER_BLOCK
569
570 #undef EL_EWALD_ANY
571 #undef EXCLUSION_FORCES
572 #undef LJ_EWALD
573
574 #undef LJ_COMB