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