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