Simple union-find structure
[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,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  */
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  *  \author Berk Hess <hess@kth.se>
45  *  \ingroup module_mdlib
46  */
47
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  */
56
57 #if GMX_PTX_ARCH < 300 && GMX_PTX_ARCH != 0
58 #error "nbnxn_cuda_kernel.cuh included with GMX_PTX_ARCH < 300 or host pass"
59 #endif
60
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
65
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
75
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
80
81 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
82 #define LJ_COMB
83 #endif
84
85 /*
86    Kernel launch parameters:
87     - #blocks   = #pair lists, blockId = pair list Id
88     - #threads  = NTHREAD_Z * c_clSize^2
89     - shmem     = see nbnxn_cuda.cu:calc_shmem_required_nonbonded()
90
91     Each thread calculates an i force-component taking one pair of i-j atoms.
92  */
93
94 /**@{*/
95 /*! \brief Compute capability dependent definition of kernel launch configuration parameters.
96  *
97  * NTHREAD_Z controls the number of j-clusters processed concurrently on NTHREAD_Z
98  * warp-pairs per block.
99  *
100  * - On CC 2.0-3.5, and >=5.0 NTHREAD_Z == 1, translating to 64 th/block with 16
101  * blocks/multiproc, is the fastest even though this setup gives low occupancy
102  * (except on 6.0).
103  * NTHREAD_Z > 1 results in excessive register spilling unless the minimum blocks
104  * per multiprocessor is reduced proportionally to get the original number of max
105  * threads in flight (and slightly lower performance).
106  * - On CC 3.7 there are enough registers to double the number of threads; using
107  * NTHREADS_Z == 2 is fastest with 16 blocks (TODO: test with RF and other kernels
108  * with low-register use).
109  *
110  * Note that the current kernel implementation only supports NTHREAD_Z > 1 with
111  * shuffle-based reduction, hence CC >= 3.0.
112  */
113
114 /* Kernel launch bounds for different compute capabilities. The value of NTHREAD_Z
115  * determines the number of threads per block and it is chosen such that
116  * 16 blocks/multiprocessor can be kept in flight.
117  * - CC 3.0,3.5, and >=5.0: NTHREAD_Z=1, (64, 16) bounds
118  * - CC 3.7:                NTHREAD_Z=2, (128, 16) bounds
119  *
120  * Note: convenience macros, need to be undef-ed at the end of the file.
121  */
122 #if GMX_PTX_ARCH == 370
123     #define NTHREAD_Z           (2)
124     #define MIN_BLOCKS_PER_MP   (16)
125 #else
126     #define NTHREAD_Z           (1)
127     #define MIN_BLOCKS_PER_MP   (16)
128 #endif /* GMX_PTX_ARCH == 370 */
129 #define THREADS_PER_BLOCK   (c_clSize*c_clSize*NTHREAD_Z)
130
131 #if GMX_PTX_ARCH >= 350
132 #if (GMX_PTX_ARCH <= 210) && (NTHREAD_Z > 1)
133     #error NTHREAD_Z > 1 will give incorrect results on CC 2.x
134 #endif
135 /**@}*/
136 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
137 #else
138 __launch_bounds__(THREADS_PER_BLOCK)
139 #endif /* GMX_PTX_ARCH >= 350 */
140 #ifdef PRUNE_NBL
141 #ifdef CALC_ENERGIES
142 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
143 #else
144 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
145 #endif /* CALC_ENERGIES */
146 #else
147 #ifdef CALC_ENERGIES
148 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
149 #else
150 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
151 #endif /* CALC_ENERGIES */
152 #endif /* PRUNE_NBL */
153 (const cu_atomdata_t atdat,
154  const cu_nbparam_t nbparam,
155  const cu_plist_t plist,
156  bool bCalcFshift)
157 #ifdef FUNCTION_DECLARATION_ONLY
158 ;     /* Only do function declaration, omit the function body. */
159 #else
160 {
161     /* convenience variables */
162     const nbnxn_sci_t *pl_sci       = plist.sci;
163 #ifndef PRUNE_NBL
164     const
165 #endif
166     nbnxn_cj4_t        *pl_cj4      = plist.cj4;
167     const nbnxn_excl_t *excl        = plist.excl;
168 #ifndef LJ_COMB
169     const int          *atom_types  = atdat.atom_types;
170     int                 ntypes      = atdat.ntypes;
171 #else
172     const float2       *lj_comb     = atdat.lj_comb;
173     float2              ljcp_i, ljcp_j;
174 #endif
175     const float4       *xq          = atdat.xq;
176     float3             *f           = atdat.f;
177     const float3       *shift_vec   = atdat.shift_vec;
178     float               rcoulomb_sq = nbparam.rcoulomb_sq;
179 #ifdef VDW_CUTOFF_CHECK
180     float               rvdw_sq     = nbparam.rvdw_sq;
181     float               vdw_in_range;
182 #endif
183 #ifdef LJ_EWALD
184     float               lje_coeff2, lje_coeff6_6;
185 #endif
186 #ifdef EL_RF
187     float two_k_rf              = nbparam.two_k_rf;
188 #endif
189 #ifdef EL_EWALD_ANA
190     float beta2                 = nbparam.ewald_beta*nbparam.ewald_beta;
191     float beta3                 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
192 #endif
193 #ifdef PRUNE_NBL
194     float rlist_sq              = nbparam.rlistOuter_sq;
195 #endif
196
197 #ifdef CALC_ENERGIES
198 #ifdef EL_EWALD_ANY
199     float  beta        = nbparam.ewald_beta;
200     float  ewald_shift = nbparam.sh_ewald;
201 #else
202     float  c_rf        = nbparam.c_rf;
203 #endif /* EL_EWALD_ANY */
204     float *e_lj        = atdat.e_lj;
205     float *e_el        = atdat.e_el;
206 #endif /* CALC_ENERGIES */
207
208     /* thread/block/warp id-s */
209     unsigned int tidxi  = threadIdx.x;
210     unsigned int tidxj  = threadIdx.y;
211     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
212 #if NTHREAD_Z == 1
213     unsigned int tidxz  = 0;
214 #else
215     unsigned int tidxz  = threadIdx.z;
216 #endif
217     unsigned int bidx   = blockIdx.x;
218     unsigned int widx   = tidx / warp_size; /* warp index */
219
220     int          sci, ci, cj,
221                  ai, aj,
222                  cij4_start, cij4_end;
223 #ifndef LJ_COMB
224     int          typei, typej;
225 #endif
226     int          i, jm, j4, wexcl_idx;
227     float        qi, qj_f,
228                  r2, inv_r, inv_r2;
229 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
230     float        inv_r6, c6, c12;
231 #endif
232 #ifdef LJ_COMB_LB
233     float        sigma, epsilon;
234 #endif
235     float        int_bit,
236                  F_invr;
237 #ifdef CALC_ENERGIES
238     float        E_lj, E_el;
239 #endif
240 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
241     float        E_lj_p;
242 #endif
243     unsigned int wexcl, imask, mask_ji;
244     float4       xqbuf;
245     float3       xi, xj, rv, f_ij, fcj_buf;
246     float3       fci_buf[c_numClPerSupercl]; /* i force buffer */
247     nbnxn_sci_t  nb_sci;
248
249     /*! i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
250     const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
251
252     /* shmem buffer for i x+q pre-loading */
253     extern __shared__  float4 xqib[];
254
255     /* shmem buffer for cj, for each warp separately */
256     int *cjs       = ((int *)(xqib + c_numClPerSupercl * c_clSize)) + tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
257     int *cjs_end   = ((int *)(xqib + c_numClPerSupercl * c_clSize)) + NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
258 #ifndef LJ_COMB
259     /* shmem buffer for i atom-type pre-loading */
260     int *atib      = cjs_end;
261 #else
262     /* shmem buffer for i-atom LJ combination rule parameters */
263     float2 *ljcpib = (float2 *)cjs_end;
264 #endif
265
266     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
267     sci         = nb_sci.sci;           /* super-cluster */
268     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
269     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
270
271     if (tidxz == 0)
272     {
273         /* Pre-load i-atom x and q into shared memory */
274         ci = sci * c_numClPerSupercl + tidxj;
275         ai = ci * c_clSize + tidxi;
276
277         float  *shiftptr = (float *)&shift_vec[nb_sci.shift];
278         xqbuf    = xq[ai] + make_float4(LDG(shiftptr), LDG(shiftptr + 1), LDG(shiftptr + 2), 0.0f);
279         xqbuf.w *= nbparam.epsfac;
280         xqib[tidxj * c_clSize + tidxi] = xqbuf;
281
282 #ifndef LJ_COMB
283         /* Pre-load the i-atom types into shared memory */
284         atib[tidxj * c_clSize + tidxi] = atom_types[ai];
285 #else
286         /* Pre-load the LJ combination parameters into shared memory */
287         ljcpib[tidxj * c_clSize + tidxi] = lj_comb[ai];
288 #endif
289     }
290     __syncthreads();
291
292     for (i = 0; i < c_numClPerSupercl; i++)
293     {
294         fci_buf[i] = make_float3(0.0f);
295     }
296
297 #ifdef LJ_EWALD
298     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
299     lje_coeff2   = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
300     lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*c_oneSixth;
301 #endif
302
303
304 #ifdef CALC_ENERGIES
305     E_lj = 0.0f;
306     E_el = 0.0f;
307
308 #ifdef EXCLUSION_FORCES /* Ewald or RF */
309     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*c_numClPerSupercl)
310     {
311         /* we have the diagonal: add the charge and LJ self interaction energy term */
312         for (i = 0; i < c_numClPerSupercl; i++)
313         {
314 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
315             qi    = xqib[i * c_clSize + tidxi].w;
316             E_el += qi*qi;
317 #endif
318
319 #ifdef LJ_EWALD
320     #if DISABLE_CUDA_TEXTURES
321             E_lj += LDG(&nbparam.nbfp[atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2]);
322     #else
323             E_lj += tex1Dfetch<float>(nbparam.nbfp_texobj, atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2);
324     #endif
325 #endif
326         }
327
328         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
329 #ifdef LJ_EWALD
330         E_lj /= c_clSize*NTHREAD_Z;
331         E_lj *= 0.5f*c_oneSixth*lje_coeff6_6;
332 #endif
333
334 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
335         /* Correct for epsfac^2 due to adding qi^2 */
336         E_el /= nbparam.epsfac*c_clSize*NTHREAD_Z;
337 #if defined EL_RF || defined EL_CUTOFF
338         E_el *= -0.5f*c_rf;
339 #else
340         E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
341 #endif
342 #endif                                  /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
343     }
344 #endif                                  /* EXCLUSION_FORCES */
345
346 #endif                                  /* CALC_ENERGIES */
347
348 #ifdef EXCLUSION_FORCES
349     const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
350 #endif
351
352     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
353     for (j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
354     {
355         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
356         imask       = pl_cj4[j4].imei[widx].imask;
357         wexcl       = excl[wexcl_idx].pair[(tidx) & (warp_size - 1)];
358
359 #ifndef PRUNE_NBL
360         if (imask)
361 #endif
362         {
363             /* Pre-load cj into shared memory on both warps separately */
364             if ((tidxj == 0 | tidxj == 4) & (tidxi < c_nbnxnGpuJgroupSize))
365             {
366                 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize/c_splitClSize] = pl_cj4[j4].cj[tidxi];
367             }
368
369             /* Unrolling this loop
370                - with pruning leads to register spilling;
371                - on Kepler and later it is much slower;
372                Tested with up to nvcc 7.5 */
373             for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
374             {
375                 if (imask & (superClInteractionMask << (jm * c_numClPerSupercl)))
376                 {
377                     mask_ji = (1U << (jm * c_numClPerSupercl));
378
379                     cj      = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize/c_splitClSize];
380                     aj      = cj * c_clSize + tidxj;
381
382                     /* load j atom data */
383                     xqbuf   = xq[aj];
384                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
385                     qj_f    = xqbuf.w;
386 #ifndef LJ_COMB
387                     typej   = atom_types[aj];
388 #else
389                     ljcp_j  = lj_comb[aj];
390 #endif
391
392                     fcj_buf = make_float3(0.0f);
393
394 #if !defined PRUNE_NBL
395 #pragma unroll 8
396 #endif
397                     for (i = 0; i < c_numClPerSupercl; i++)
398                     {
399                         if (imask & mask_ji)
400                         {
401                             ci      = sci * c_numClPerSupercl + i; /* i cluster index */
402
403                             /* all threads load an atom from i cluster ci into shmem! */
404                             xqbuf   = xqib[i * c_clSize + tidxi];
405                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
406
407                             /* distance between i and j atoms */
408                             rv      = xi - xj;
409                             r2      = norm2(rv);
410
411 #ifdef PRUNE_NBL
412                             /* If _none_ of the atoms pairs are in cutoff range,
413                                the bit corresponding to the current
414                                cluster-pair in imask gets set to 0. */
415                             if (!__any(r2 < rlist_sq))
416                             {
417                                 imask &= ~mask_ji;
418                             }
419 #endif
420
421                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
422
423                             /* cutoff & exclusion check */
424 #ifdef EXCLUSION_FORCES
425                             if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
426 #else
427                             if ((r2 < rcoulomb_sq) * int_bit)
428 #endif
429                             {
430                                 /* load the rest of the i-atom parameters */
431                                 qi      = xqbuf.w;
432
433 #ifndef LJ_COMB
434                                 /* LJ 6*C6 and 12*C12 */
435                                 typei   = atib[i * c_clSize + tidxi];
436                                 fetch_nbfp_c6_c12(c6, c12, nbparam, ntypes * typei + typej);
437 #else
438                                 ljcp_i  = ljcpib[i * c_clSize + tidxi];
439 #ifdef LJ_COMB_GEOM
440                                 c6      = ljcp_i.x * ljcp_j.x;
441                                 c12     = ljcp_i.y * ljcp_j.y;
442 #else
443                                 /* LJ 2^(1/6)*sigma and 12*epsilon */
444                                 sigma   = ljcp_i.x + ljcp_j.x;
445                                 epsilon = ljcp_i.y * ljcp_j.y;
446 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
447                                 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
448 #endif
449 #endif                          /* LJ_COMB_GEOM */
450 #endif                          /* LJ_COMB */
451
452                                 // Ensure distance do not become so small that r^-12 overflows
453                                 r2      = max(r2, NBNXN_MIN_RSQ);
454
455                                 inv_r   = rsqrt(r2);
456                                 inv_r2  = inv_r * inv_r;
457 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
458                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
459 #ifdef EXCLUSION_FORCES
460                                 /* We could mask inv_r2, but with Ewald
461                                  * masking both inv_r6 and F_invr is faster */
462                                 inv_r6  *= int_bit;
463 #endif                          /* EXCLUSION_FORCES */
464
465                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
466 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
467                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*c_oneTwelveth -
468                                                      c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*c_oneSixth);
469 #endif
470 #else                           /* !LJ_COMB_LB || CALC_ENERGIES */
471                                 float sig_r  = sigma*inv_r;
472                                 float sig_r2 = sig_r*sig_r;
473                                 float sig_r6 = sig_r2*sig_r2*sig_r2;
474 #ifdef EXCLUSION_FORCES
475                                 sig_r6 *= int_bit;
476 #endif                          /* EXCLUSION_FORCES */
477
478                                 F_invr  = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
479 #endif                          /* !LJ_COMB_LB || CALC_ENERGIES */
480
481 #ifdef LJ_FORCE_SWITCH
482 #ifdef CALC_ENERGIES
483                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
484 #else
485                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
486 #endif /* CALC_ENERGIES */
487 #endif /* LJ_FORCE_SWITCH */
488
489
490 #ifdef LJ_EWALD
491 #ifdef LJ_EWALD_COMB_GEOM
492 #ifdef CALC_ENERGIES
493                                 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);
494 #else
495                                 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
496 #endif                          /* CALC_ENERGIES */
497 #elif defined LJ_EWALD_COMB_LB
498                                 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
499 #ifdef CALC_ENERGIES
500                                                                int_bit, &F_invr, &E_lj_p
501 #else
502                                                                0, &F_invr, NULL
503 #endif /* CALC_ENERGIES */
504                                                                );
505 #endif /* LJ_EWALD_COMB_GEOM */
506 #endif /* LJ_EWALD */
507
508 #ifdef LJ_POT_SWITCH
509 #ifdef CALC_ENERGIES
510                                 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
511 #else
512                                 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
513 #endif /* CALC_ENERGIES */
514 #endif /* LJ_POT_SWITCH */
515
516 #ifdef VDW_CUTOFF_CHECK
517                                 /* Separate VDW cut-off check to enable twin-range cut-offs
518                                  * (rvdw < rcoulomb <= rlist)
519                                  */
520                                 vdw_in_range  = (r2 < rvdw_sq) ? 1.0f : 0.0f;
521                                 F_invr       *= vdw_in_range;
522 #ifdef CALC_ENERGIES
523                                 E_lj_p       *= vdw_in_range;
524 #endif
525 #endif                          /* VDW_CUTOFF_CHECK */
526
527 #ifdef CALC_ENERGIES
528                                 E_lj    += E_lj_p;
529 #endif
530
531
532 #ifdef EL_CUTOFF
533 #ifdef EXCLUSION_FORCES
534                                 F_invr  += qi * qj_f * int_bit * inv_r2 * inv_r;
535 #else
536                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
537 #endif
538 #endif
539 #ifdef EL_RF
540                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
541 #endif
542 #if defined EL_EWALD_ANA
543                                 F_invr  += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
544 #elif defined EL_EWALD_TAB
545                                 F_invr  += qi * qj_f * (int_bit*inv_r2 -
546                                                         interpolate_coulomb_force_r(nbparam, r2 * inv_r)) * inv_r;
547 #endif                          /* EL_EWALD_ANA/TAB */
548
549 #ifdef CALC_ENERGIES
550 #ifdef EL_CUTOFF
551                                 E_el    += qi * qj_f * (int_bit*inv_r - c_rf);
552 #endif
553 #ifdef EL_RF
554                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
555 #endif
556 #ifdef EL_EWALD_ANY
557                                 /* 1.0f - erff is faster than erfcf */
558                                 E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
559 #endif                          /* EL_EWALD_ANY */
560 #endif
561                                 f_ij    = rv * F_invr;
562
563                                 /* accumulate j forces in registers */
564                                 fcj_buf -= f_ij;
565
566                                 /* accumulate i forces in registers */
567                                 fci_buf[i] += f_ij;
568                             }
569                         }
570
571                         /* shift the mask bit by 1 */
572                         mask_ji += mask_ji;
573                     }
574
575                     /* reduce j forces */
576                     reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
577                 }
578             }
579 #ifdef PRUNE_NBL
580             /* Update the imask with the new one which does not contain the
581                out of range clusters anymore. */
582             pl_cj4[j4].imei[widx].imask = imask;
583 #endif
584         }
585     }
586
587     /* skip central shifts when summing shift forces */
588     if (nb_sci.shift == CENTRAL)
589     {
590         bCalcFshift = false;
591     }
592
593     float fshift_buf = 0.0f;
594
595     /* reduce i forces */
596     for (i = 0; i < c_numClPerSupercl; i++)
597     {
598         ai  = (sci * c_numClPerSupercl + i) * c_clSize + tidxi;
599         reduce_force_i_warp_shfl(fci_buf[i], f,
600                                  &fshift_buf, bCalcFshift,
601                                  tidxj, ai);
602     }
603
604     /* add up local shift forces into global mem, tidxj indexes x,y,z */
605     if (bCalcFshift && (tidxj & 3) < 3)
606     {
607         atomicAdd(&(atdat.fshift[nb_sci.shift].x) + (tidxj & 3), fshift_buf);
608     }
609
610 #ifdef CALC_ENERGIES
611     /* reduce the energies over warps and store into global memory */
612     reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
613 #endif
614 }
615 #endif /* FUNCTION_DECLARATION_ONLY */
616
617 #undef NTHREAD_Z
618 #undef MIN_BLOCKS_PER_MP
619 #undef THREADS_PER_BLOCK
620
621 #undef EL_EWALD_ANY
622 #undef EXCLUSION_FORCES
623 #undef LJ_EWALD
624
625 #undef LJ_COMB