Update some nbnxm kernel constants to constexpr
[alexxy/gromacs.git] / src / gromacs / nbnxm / opencl / nbnxm_ocl_kernel.clh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012-2018, The GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *  \brief OpenCL non-bonded kernel.
39  *
40  *  OpenCL 1.2 support is expected.
41  *
42  *  \author Anca Hamuraru <anca@streamcomputing.eu>
43  *  \author Szilárd Páll <pall.szilard@gmail.com>
44  *  \ingroup module_nbnxm
45  */
46
47 /* Currently we enable CJ prefetch for AMD/NVIDIA and disable it for the "nowarp" kernel
48  * Note that this should precede the kernel_utils include.
49  */
50 #include "nbnxm_ocl_kernel_utils.clh"
51
52 /////////////////////////////////////////////////////////////////////////////////////////////////
53
54 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
55 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
56 #    define EL_EWALD_ANY
57 #endif
58
59 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD \
60         || (defined EL_CUTOFF && defined CALC_ENERGIES)
61 /* Macro to control the calculation of exclusion forces in the kernel
62  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
63  * energy terms.
64  *
65  * Note: convenience macro, needs to be undef-ed at the end of the file.
66  */
67 #    define EXCLUSION_FORCES
68 #endif
69
70 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
71 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
72 #    define LJ_EWALD
73 #endif
74
75 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
76 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
77 #    define LJ_COMB
78 #endif
79
80 /*
81    Kernel launch parameters:
82     - #blocks   = #pair lists, blockId = pair list Id
83     - #threads  = CL_SIZE^2
84     - shmem     = CL_SIZE^2 * sizeof(float)
85
86     Each thread calculates an i force-component taking one pair of i-j atoms.
87
88    TODO: implement 128 threads/wavefront by porting over the NTHREAD_Z/j4 loop
89    "horizontal splitting" over threads.
90  */
91
92 /* NOTE:
93    NB_KERNEL_FUNC_NAME differs from the CUDA equivalent as it is not a variadic macro due to OpenCL
94    not having a support for them, this version only takes exactly 2 arguments. Thus if more strings
95    need to be appended a new macro must be written or it must be directly appended here.
96  */
97 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
98 #ifdef cl_intel_required_subgroup_size
99 __attribute__((intel_reqd_sub_group_size(SUBGROUP_SIZE)))
100 #endif
101 #ifdef PRUNE_NBL
102 #    ifdef CALC_ENERGIES
103 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
104 #    else
105 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
106 #    endif
107 #else
108 #    ifdef CALC_ENERGIES
109 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
110 #    else
111 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
112 #    endif
113 #endif
114         (
115 #ifndef LJ_COMB
116                 int ntypes, /* IN  */
117 #endif
118                 cl_nbparam_params_t nbparam_params,       /* IN  */
119                 const __global float4* restrict xq,       /* IN  */
120                 __global float* restrict f,               /* OUT stores float3 values */
121                 __global float* restrict gmx_unused e_lj, /* OUT */
122                 __global float* restrict gmx_unused e_el, /* OUT */
123                 __global float* restrict fshift,          /* OUT stores float3 values */
124 #ifdef LJ_COMB
125                 const __global float2* restrict lj_comb, /* IN stores float2 values */
126 #else
127                 const __global int* restrict atom_types, /* IN  */
128 #endif
129                 const __global float* restrict shift_vec,             /* IN stores float3 values */
130                 __constant const float* gmx_unused nbfp_climg2d,      /* IN  */
131                 __constant const float* gmx_unused nbfp_comb_climg2d, /* IN  */
132                 __constant const float* gmx_unused coulomb_tab_climg2d, /* IN  */
133                 const __global nbnxn_sci_t* pl_sci,                     /* IN  */
134 #ifndef PRUNE_NBL
135                 const
136 #endif
137                 __global nbnxn_cj4_t* pl_cj4,             /* OUT / IN */
138                 const __global nbnxn_excl_t* excl,        /* IN  */
139                 int                          bCalcFshift, /* IN  */
140                 __local float4* xqib                      /* Pointer to dyn alloc'ed shmem */
141         )
142 {
143     /* convenience variables */
144     const cl_nbparam_params_t* const nbparam = &nbparam_params;
145
146     const float rcoulomb_sq = nbparam->rcoulomb_sq;
147 #ifdef VDW_CUTOFF_CHECK
148     const float rvdw_sq = nbparam_params.rvdw_sq;
149 #endif
150 #ifdef EL_RF
151     const float two_k_rf = nbparam->two_k_rf;
152 #endif
153 #ifdef EL_EWALD_TAB
154     const float coulomb_tab_scale = nbparam->coulomb_tab_scale;
155 #endif
156 #ifdef EL_EWALD_ANA
157     const float beta2 = nbparam->ewald_beta * nbparam->ewald_beta;
158     const float beta3 = nbparam->ewald_beta * nbparam->ewald_beta * nbparam->ewald_beta;
159 #endif
160 #ifdef PRUNE_NBL
161     const float rlist_sq = nbparam->rlistOuter_sq;
162 #endif
163
164 #ifdef CALC_ENERGIES
165 #    ifdef EL_EWALD_ANY
166     const float beta        = nbparam->ewald_beta;
167     const float ewald_shift = nbparam->sh_ewald;
168 #    else
169     const float gmx_unused c_rf = nbparam->c_rf;
170 #    endif /* EL_EWALD_ANY */
171 #endif     /* CALC_ENERGIES */
172
173     /* thread/block/warp id-s */
174     const int tidxi = get_local_id(0);
175     const int tidxj = get_local_id(1);
176     const int tidx  = (int)(get_local_id(1) * get_local_size(0) + get_local_id(0));
177     const int bidx  = get_group_id(0);
178     const int widx  = tidx / WARP_SIZE; /* warp index */
179
180     /*! i-cluster interaction mask for a super-cluster with all c_nbnxnGpuNumClusterPerSupercluster=8 bits set */
181     const unsigned superClInteractionMask = ((1U << c_nbnxnGpuNumClusterPerSupercluster) - 1U);
182
183 #define LOCAL_OFFSET (xqib + c_nbnxnGpuNumClusterPerSupercluster * CL_SIZE)
184     CjType cjs = 0;
185 #if USE_CJ_PREFETCH
186     /* shmem buffer for cj, for both warps separately */
187     cjs = (__local int*)(LOCAL_OFFSET);
188 #    undef LOCAL_OFFSET
189 #    define LOCAL_OFFSET (cjs + 2 * c_nbnxnGpuJgroupSize)
190 #endif // USE_CJ_PREFETCH
191
192 #ifdef IATYPE_SHMEM
193 #    ifndef LJ_COMB
194     /* shmem buffer for i atom-type pre-loading */
195     __local int* atib = (__local int*)(LOCAL_OFFSET); //NOLINT(google-readability-casting)
196 #        undef LOCAL_OFFSET
197 #        define LOCAL_OFFSET (atib + c_nbnxnGpuNumClusterPerSupercluster * CL_SIZE)
198 #    else
199     __local float2* ljcpib      = (__local float2*)(LOCAL_OFFSET);
200 #        undef LOCAL_OFFSET
201 #        define LOCAL_OFFSET (ljcpib + c_nbnxnGpuNumClusterPerSupercluster * CL_SIZE)
202 #    endif
203 #endif
204
205 #if !REDUCE_SHUFFLE
206     /* shmem j force buffer */
207     __local float* f_buf = (__local float*)(LOCAL_OFFSET);
208 #    undef LOCAL_OFFSET
209 #    define LOCAL_OFFSET (f_buf + CL_SIZE * CL_SIZE * 3)
210 #else
211     __local float* f_buf             = 0;
212 #endif
213 #if !USE_SUBGROUP_ANY
214     /* Local buffer used to implement __any warp vote function from CUDA.
215        volatile is used to avoid compiler optimizations for AMD builds. */
216     //NOLINTNEXTLINE(google-readability-casting)
217     volatile __local int* warp_any = (__local int*)(LOCAL_OFFSET);
218 #else
219     __local int gmx_unused* warp_any = 0;
220 #endif
221 #undef LOCAL_OFFSET
222
223     const nbnxn_sci_t nb_sci     = pl_sci[bidx]; /* my i super-cluster's index = current bidx */
224     const int         sci        = nb_sci.sci;   /* super-cluster */
225     const int         cij4_start = nb_sci.cj4_ind_start; /* first ...*/
226     const int         cij4_end   = nb_sci.cj4_ind_end;   /* and last index of j clusters */
227
228     for (int i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i += CL_SIZE)
229     {
230         /* Pre-load i-atom x and q into shared memory */
231         const int ci = sci * c_nbnxnGpuNumClusterPerSupercluster + tidxj + i;
232         const int ai = ci * CL_SIZE + tidxi;
233
234         float4 xqbuf = xq[ai]
235                        + (float4)(shift_vec[3 * nb_sci.shift], shift_vec[3 * nb_sci.shift + 1],
236                                   shift_vec[3 * nb_sci.shift + 2], 0.0F);
237         xqbuf.w *= nbparam->epsfac;
238         xqib[(tidxj + i) * CL_SIZE + tidxi] = xqbuf;
239 #ifdef IATYPE_SHMEM
240 #    ifndef LJ_COMB
241         /* Pre-load the i-atom types into shared memory */
242         atib[(tidxj + i) * CL_SIZE + tidxi] = atom_types[ai];
243 #    else
244         ljcpib[(tidxj + i) * CL_SIZE + tidxi] = lj_comb[ai];
245 #    endif
246 #endif
247     }
248 #if !USE_SUBGROUP_ANY
249     /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
250     if (tidx == 0 || tidx == WARP_SIZE)
251     {
252         warp_any[widx] = 0;
253     }
254 #endif
255     barrier(CLK_LOCAL_MEM_FENCE);
256
257     float3 fci_buf[c_nbnxnGpuNumClusterPerSupercluster]; /* i force buffer */
258     for (int ci_offset = 0; ci_offset < c_nbnxnGpuNumClusterPerSupercluster; ci_offset++)
259     {
260         fci_buf[ci_offset] = (float3)(0.0F);
261     }
262
263 #ifdef LJ_EWALD
264     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
265     const float lje_coeff2   = nbparam->ewaldcoeff_lj * nbparam->ewaldcoeff_lj;
266     const float lje_coeff6_6 = lje_coeff2 * lje_coeff2 * lje_coeff2 * ONE_SIXTH_F;
267 #endif /* LJ_EWALD */
268
269
270 #ifdef CALC_ENERGIES
271     float E_lj = 0.0F;
272     float E_el = 0.0F;
273
274 #    if defined EXCLUSION_FORCES /* Ewald or RF */
275     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci * c_nbnxnGpuNumClusterPerSupercluster)
276     {
277         /* we have the diagonal: add the charge and LJ self interaction energy term */
278         for (int i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
279         {
280 #        if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
281             const float qi = xqib[i * CL_SIZE + tidxi].w;
282             E_el += qi * qi;
283 #        endif
284 #        if defined LJ_EWALD
285             E_lj += nbfp_climg2d[atom_types[(sci * c_nbnxnGpuNumClusterPerSupercluster + i) * CL_SIZE + tidxi]
286                                  * (ntypes + 1) * 2];
287 #        endif /* LJ_EWALD */
288         }
289
290         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
291 #        ifdef LJ_EWALD
292         E_lj /= CL_SIZE;
293         E_lj *= HALF_F * ONE_SIXTH_F * lje_coeff6_6;
294 #        endif /* LJ_EWALD */
295
296 #        if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
297         /* Correct for epsfac^2 due to adding qi^2 */
298         E_el /= nbparam->epsfac * CL_SIZE;
299 #            if defined EL_RF || defined EL_CUTOFF
300         E_el *= -HALF_F * c_rf;
301 #            else
302         E_el *= -beta * M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
303 #            endif
304 #        endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
305     }
306 #    endif /* EXCLUSION_FORCES */
307
308 #endif /* CALC_ENERGIES */
309
310 #ifdef EXCLUSION_FORCES
311     const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
312 #endif
313
314     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
315     for (int j4 = cij4_start; j4 < cij4_end; j4++)
316     {
317         const int          wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
318         unsigned int       imask     = pl_cj4[j4].imei[widx].imask;
319         const unsigned int wexcl     = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
320
321         preloadCj4(&cjs, pl_cj4[j4].cj, tidxi, tidxj, imask != 0U);
322
323 #ifndef PRUNE_NBL
324         if (imask)
325 #endif
326         {
327             /* Unrolling this loop improves performance without pruning but
328              * with pruning it leads to slowdown.
329              *
330              * Tested with driver 1800.5
331              *
332              * TODO: check loop unrolling with NVIDIA OpenCL
333              */
334 #if !defined PRUNE_NBL && !defined _NVIDIA_SOURCE_
335 #    pragma unroll 4
336 #endif
337             for (int jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
338             {
339                 if (imask & (superClInteractionMask << (jm * c_nbnxnGpuNumClusterPerSupercluster)))
340                 {
341                     unsigned int mask_ji = (1U << (jm * c_nbnxnGpuNumClusterPerSupercluster));
342
343                     const int cj = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
344                     const int aj = cj * CL_SIZE + tidxj;
345
346                     /* load j atom data */
347                     const float4 xjqbuf = xq[aj];
348                     const float3 xj     = (float3)(xjqbuf.xyz);
349                     const float  qj_f   = xjqbuf.w;
350 #ifndef LJ_COMB
351                     const int typej = atom_types[aj];
352 #else
353                     const float2 ljcp_j = lj_comb[aj];
354 #endif
355
356                     float3 fcj_buf = (float3)(0.0F);
357
358 #if !defined PRUNE_NBL
359 #    pragma unroll 8
360 #endif
361                     for (int i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
362                     {
363                         if (imask & mask_ji)
364                         {
365                             const int gmx_unused ci = sci * c_nbnxnGpuNumClusterPerSupercluster + i; /* i cluster index */
366
367                             /* all threads load an atom from i cluster ci into shmem! */
368                             const float4 xiqbuf = xqib[i * CL_SIZE + tidxi];
369                             const float3 xi     = (float3)(xiqbuf.xyz);
370
371                             /* distance between i and j atoms */
372                             const float3 rv = xi - xj;
373                             float        r2 = norm2(rv);
374
375 #ifdef PRUNE_NBL
376                             if (!gmx_sub_group_any(warp_any, widx, r2 < rlist_sq))
377                             {
378                                 imask &= ~mask_ji;
379                             }
380 #endif
381
382                             const float int_bit = (wexcl & mask_ji) ? 1.0F : 0.0F;
383
384                             /* cutoff & exclusion check */
385 #ifdef EXCLUSION_FORCES
386                             if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
387 #else
388                             if ((float)(r2 < rcoulomb_sq) * int_bit != 0.0F)
389 #endif
390                             {
391                                 /* load the rest of the i-atom parameters */
392                                 const float qi = xiqbuf.w;
393 #ifdef IATYPE_SHMEM
394 #    ifndef LJ_COMB
395                                 const int typei = atib[i * CL_SIZE + tidxi];
396 #    else
397                                 const float2 ljcp_i = ljcpib[i * CL_SIZE + tidxi];
398 #    endif
399 #else /* IATYPE_SHMEM */
400                                 const int   ai     = ci * CL_SIZE + tidxi; /* i atom index */
401
402 #    ifndef LJ_COMB
403                                 const int   typei  = atom_types[ai];
404 #    else
405                                 const float2 ljcp_i = lj_comb[ai];
406 #    endif
407 #endif                          /* IATYPE_SHMEM */
408                                 /* LJ 6*C6 and 12*C12 */
409 #ifndef LJ_COMB
410                                 const float c6  = nbfp_climg2d[2 * (ntypes * typei + typej)];
411                                 const float c12 = nbfp_climg2d[2 * (ntypes * typei + typej) + 1];
412 #else /* LJ_COMB */
413 #    ifdef LJ_COMB_GEOM
414                                 const float c6     = ljcp_i.x * ljcp_j.x;
415                                 const float c12    = ljcp_i.y * ljcp_j.y;
416 #    else
417                                 /* LJ 2^(1/6)*sigma and 12*epsilon */
418                                 const float sigma   = ljcp_i.x + ljcp_j.x;
419                                 const float epsilon = ljcp_i.y * ljcp_j.y;
420 #        if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
421                                 float       c6, c12;
422                                 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
423 #        endif
424 #    endif /* LJ_COMB_GEOM */
425 #endif     /* LJ_COMB */
426
427                                 // Ensure distance do not become so small that r^-12 overflows.
428                                 // Cast to float to ensure the correct built-in max() function
429                                 // is called.
430                                 r2 = max(r2, (float)c_nbnxnMinDistanceSquared);
431
432                                 const float inv_r  = rsqrt(r2);
433                                 const float inv_r2 = inv_r * inv_r;
434 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
435                                 float inv_r6 = inv_r2 * inv_r2 * inv_r2;
436 #    if defined EXCLUSION_FORCES
437                                 /* We could mask inv_r2, but with Ewald
438                                  * masking both inv_r6 and F_invr is faster */
439                                 inv_r6 *= int_bit;
440 #    endif /* EXCLUSION_FORCES */
441
442                                 float F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
443 #    if defined CALC_ENERGIES || defined LJ_POT_SWITCH
444                                 float E_lj_p =
445                                         int_bit
446                                         * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot) * ONE_TWELVETH_F
447                                            - c6 * (inv_r6 + nbparam->dispersion_shift.cpot) * ONE_SIXTH_F);
448
449 #    endif
450 #else /* ! LJ_COMB_LB || CALC_ENERGIES */
451                                 const float sig_r  = sigma * inv_r;
452                                 const float sig_r2 = sig_r * sig_r;
453                                 float       sig_r6 = sig_r2 * sig_r2 * sig_r2;
454 #    if defined EXCLUSION_FORCES
455                                 sig_r6 *= int_bit;
456 #    endif /* EXCLUSION_FORCES */
457
458                                 float F_invr = epsilon * sig_r6 * (sig_r6 - 1.0F) * inv_r2;
459 #endif     /* ! LJ_COMB_LB || CALC_ENERGIES */
460
461
462 #ifdef LJ_FORCE_SWITCH
463 #    ifdef CALC_ENERGIES
464                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
465 #    else
466                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
467 #    endif /* CALC_ENERGIES */
468 #endif     /* LJ_FORCE_SWITCH */
469
470
471 #ifdef LJ_EWALD
472 #    ifdef LJ_EWALD_COMB_GEOM
473 #        ifdef CALC_ENERGIES
474                                 calculate_lj_ewald_comb_geom_F_E(
475                                         nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2,
476                                         lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
477 #        else
478                                 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2,
479                                                                lje_coeff2, lje_coeff6_6, &F_invr);
480 #        endif /* CALC_ENERGIES */
481 #    elif defined LJ_EWALD_COMB_LB
482                                 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej,
483                                                                r2, inv_r2, lje_coeff2, lje_coeff6_6,
484 #        ifdef CALC_ENERGIES
485                                                                int_bit, true, &F_invr, &E_lj_p
486 #        else
487                                                                0, false, &F_invr, 0
488 #        endif /* CALC_ENERGIES */
489                                 );
490 #    endif     /* LJ_EWALD_COMB_GEOM */
491 #endif         /* LJ_EWALD */
492
493 #ifdef LJ_POT_SWITCH
494 #    ifdef CALC_ENERGIES
495                                 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
496 #    else
497                                 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
498 #    endif /* CALC_ENERGIES */
499 #endif     /* LJ_POT_SWITCH */
500
501 #ifdef VDW_CUTOFF_CHECK
502                                 /* Separate VDW cut-off check to enable twin-range cut-offs
503                                  * (rvdw < rcoulomb <= rlist)
504                                  */
505                                 const float vdw_in_range = (r2 < rvdw_sq) ? 1.0F : 0.0F;
506                                 F_invr *= vdw_in_range;
507 #    ifdef CALC_ENERGIES
508                                 E_lj_p *= vdw_in_range;
509 #    endif
510 #endif /* VDW_CUTOFF_CHECK */
511
512 #ifdef CALC_ENERGIES
513                                 E_lj += E_lj_p;
514
515 #endif
516
517
518 #ifdef EL_CUTOFF
519 #    ifdef EXCLUSION_FORCES
520                                 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
521 #    else
522                                 F_invr += qi * qj_f * inv_r2 * inv_r;
523 #    endif
524 #endif
525 #ifdef EL_RF
526                                 F_invr += qi * qj_f * (int_bit * inv_r2 * inv_r - two_k_rf);
527 #endif
528 #if defined EL_EWALD_ANA
529                                 F_invr += qi * qj_f
530                                           * (int_bit * inv_r2 * inv_r + pmecorrF(beta2 * r2) * beta3);
531 #elif defined EL_EWALD_TAB
532                                 F_invr += qi * qj_f
533                                           * (int_bit * inv_r2
534                                              - interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r,
535                                                                            coulomb_tab_scale))
536                                           * inv_r;
537 #endif /* EL_EWALD_ANA/TAB */
538
539 #ifdef CALC_ENERGIES
540 #    ifdef EL_CUTOFF
541                                 E_el += qi * qj_f * (int_bit * inv_r - c_rf);
542 #    endif
543 #    ifdef EL_RF
544                                 E_el += qi * qj_f * (int_bit * inv_r + HALF_F * two_k_rf * r2 - c_rf);
545 #    endif
546 #    ifdef EL_EWALD_ANY
547                                 /* 1.0F - erff is faster than erfcf */
548                                 E_el += qi * qj_f
549                                         * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
550 #    endif /* EL_EWALD_ANY */
551 #endif
552                                 const float3 f_ij = rv * F_invr;
553
554                                 /* accumulate j forces in registers */
555                                 fcj_buf -= f_ij;
556
557                                 /* accumulate i forces in registers */
558                                 fci_buf[i] += f_ij;
559                             }
560                         }
561
562                         /* shift the mask bit by 1 */
563                         mask_ji += mask_ji;
564                     }
565
566                     /* reduce j forces */
567                     reduce_force_j(f_buf, fcj_buf, f, tidxi, tidxj, aj);
568                 }
569             }
570 #ifdef PRUNE_NBL
571             /* Update the imask with the new one which does not contain the
572                out of range clusters anymore. */
573
574             pl_cj4[j4].imei[widx].imask = imask;
575 #endif
576         }
577     }
578
579     /* skip central shifts when summing shift forces */
580     if (nb_sci.shift == CENTRAL)
581     {
582         bCalcFshift = 0;
583     }
584     /* reduce i forces */
585     reduce_force_i_and_shift(f_buf, fci_buf, f, bCalcFshift != 0, tidxi, tidxj, sci, nb_sci.shift, fshift);
586
587 #ifdef CALC_ENERGIES
588     reduce_energy(f_buf, E_lj, E_el, e_lj, e_el, tidx);
589 #endif
590 }
591
592 #undef EL_EWALD_ANY
593 #undef EXCLUSION_FORCES
594 #undef LJ_EWALD
595
596 #undef LJ_COMB
597
598 #undef USE_CJ_PREFETCH