Use subgroup for warp_any and CJ4 prefetch
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_ocl / nbnxn_ocl_kernel.clh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2016,2017,2018, 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 OpenCL non-bonded kernel.
38  *
39  *  OpenCL 1.2 support is expected.
40  *
41  *  \author Anca Hamuraru <anca@streamcomputing.eu>
42  *  \author Szilárd Páll <pall.szilard@gmail.com>
43  *  \ingroup module_mdlib
44  */
45
46 /* Currently we enable CJ prefetch for AMD/NVIDIA and disable it for the "nowarp" kernel
47  * Note that this should precede the kernel_utils include.
48  */
49 #include "nbnxn_ocl_kernel_utils.clh"
50
51 /////////////////////////////////////////////////////////////////////////////////////////////////
52
53 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
54 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
55 #define EL_EWALD_ANY
56 #endif
57
58 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
59 /* Macro to control the calculation of exclusion forces in the kernel
60  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
61  * energy terms.
62  *
63  * Note: convenience macro, needs to be undef-ed at the end of the file.
64  */
65 #define EXCLUSION_FORCES
66 #endif
67
68 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
69 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
70 #define LJ_EWALD
71 #endif
72
73 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
74 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
75 #define LJ_COMB
76 #endif
77
78 /*
79    Kernel launch parameters:
80     - #blocks   = #pair lists, blockId = pair list Id
81     - #threads  = CL_SIZE^2
82     - shmem     = CL_SIZE^2 * sizeof(float)
83
84     Each thread calculates an i force-component taking one pair of i-j atoms.
85
86    TODO: implement 128 threads/wavefront by porting over the NTHREAD_Z/j4 loop
87    "horizontal splitting" over threads.
88  */
89
90 /* NOTE:
91    NB_KERNEL_FUNC_NAME differs from the CUDA equivalent as it is not a variadic macro due to OpenCL not having a support for them, this version only takes exactly 2 arguments.
92    Thus if more strings need to be appended a new macro must be written or it must be directly appended here.
93  */
94 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
95 #ifdef cl_intel_required_subgroup_size
96 __attribute__((intel_reqd_sub_group_size(SUBGROUP_SIZE)))
97 #endif
98 #ifdef PRUNE_NBL
99     #ifdef CALC_ENERGIES
100 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
101     #else
102 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
103     #endif
104 #else
105     #ifdef CALC_ENERGIES
106 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
107     #else
108 __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
109     #endif
110 #endif
111 (
112 #ifndef LJ_COMB
113     int ntypes,                               /* IN  */
114 #endif
115     cl_nbparam_params_t nbparam_params,       /* IN  */
116     const __global float4 *restrict xq,       /* IN  */
117     __global float *restrict f,               /* OUT stores float3 values */
118     __global float *restrict e_lj,            /* OUT */
119     __global float *restrict e_el,            /* OUT */
120     __global float *restrict fshift,          /* OUT stores float3 values */
121 #ifdef LJ_COMB
122     const __global float2 *restrict lj_comb,  /* IN stores float2 values */
123 #else
124     const __global int *restrict atom_types,  /* IN  */
125 #endif
126     const __global float *restrict shift_vec, /* IN stores float3 values */
127     __constant float* nbfp_climg2d,           /* IN  */
128     __constant float* nbfp_comb_climg2d,      /* IN  */
129     __constant float* coulomb_tab_climg2d,    /* IN  */
130     const __global nbnxn_sci_t* pl_sci,       /* IN  */
131 #ifndef PRUNE_NBL
132     const
133 #endif
134     __global nbnxn_cj4_t* pl_cj4,             /* OUT / IN */
135     const __global nbnxn_excl_t* excl,        /* IN  */
136     int bCalcFshift,                          /* IN  */
137     __local  float4   *xqib                   /* Pointer to dyn alloc'ed shmem */
138 )
139 {
140     /* convenience variables */
141     cl_nbparam_params_t *nbparam = &nbparam_params;
142
143     float                rcoulomb_sq = nbparam->rcoulomb_sq;
144 #ifdef LJ_COMB
145     float2               ljcp_i, ljcp_j;
146 #endif
147 #ifdef VDW_CUTOFF_CHECK
148     float               rvdw_sq     = nbparam_params.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->rlistOuter_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 #endif /* CALC_ENERGIES */
176
177     /* thread/block/warp id-s */
178     unsigned int tidxi  = get_local_id(0);
179     unsigned int tidxj  = get_local_id(1);
180     unsigned int tidx   = get_local_id(1) * get_local_size(0) + get_local_id(0);
181     unsigned int bidx   = get_group_id(0);
182     unsigned int widx   = tidx / WARP_SIZE; /* warp index */
183     int          sci, ci, cj, ci_offset,
184                  ai, aj,
185                  cij4_start, cij4_end;
186 #ifndef LJ_COMB
187     int          typei, typej;
188 #endif
189     int          i, jm, j4, wexcl_idx;
190     float        qi, qj_f,
191                  r2, inv_r, inv_r2;
192 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
193     float        inv_r6, c6, c12;
194 #endif
195 #if defined LJ_COMB_LB
196     float        sigma, epsilon;
197 #endif
198     float        int_bit,
199                  F_invr;
200
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     /*! i-cluster interaction mask for a super-cluster with all NCL_PER_SUPERCL=8 bits set */
214     const unsigned superClInteractionMask = ((1U << NCL_PER_SUPERCL) - 1U);
215
216 #define LOCAL_OFFSET (xqib + NCL_PER_SUPERCL * CL_SIZE)
217     CjType cjs;
218 #if USE_CJ_PREFETCH
219     /* shmem buffer for cj, for both warps separately */
220     cjs = (__local int *)(LOCAL_OFFSET);
221     #undef LOCAL_OFFSET
222     #define LOCAL_OFFSET cjs + 2 * NBNXN_GPU_JGROUP_SIZE
223 #endif //USE_CJ_PREFETCH
224
225 #ifdef IATYPE_SHMEM
226 #ifndef LJ_COMB
227     /* shmem buffer for i atom-type pre-loading */
228     __local int *atib      = (__local int *)(LOCAL_OFFSET);
229     #undef LOCAL_OFFSET
230     #define LOCAL_OFFSET atib + NCL_PER_SUPERCL * CL_SIZE
231 #else
232     __local float2 *ljcpib = (__local float2 *)(LOCAL_OFFSET);
233     #undef LOCAL_OFFSET
234     #define LOCAL_OFFSET ljcpib + NCL_PER_SUPERCL * CL_SIZE
235 #endif
236 #endif
237
238 #if !REDUCE_SHUFFLE
239     /* shmem j force buffer */
240     __local float *f_buf   = (__local float *)(LOCAL_OFFSET);
241     #undef LOCAL_OFFSET
242     #define LOCAL_OFFSET (f_buf + CL_SIZE * CL_SIZE * 3)
243 #else
244     __local float *f_buf  = 0;
245 #endif
246 #if !USE_SUBGROUP_ANY
247     /* Local buffer used to implement __any warp vote function from CUDA.
248        volatile is used to avoid compiler optimizations for AMD builds. */
249     volatile __local uint *warp_any = (__local uint*)(LOCAL_OFFSET);
250 #else
251     __local uint          *warp_any = 0;
252 #endif
253 #undef LOCAL_OFFSET
254
255     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
256     sci         = nb_sci.sci;           /* super-cluster */
257     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
258     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
259
260     for (i = 0; i < NCL_PER_SUPERCL; i += CL_SIZE)
261     {
262         /* Pre-load i-atom x and q into shared memory */
263         ci = sci * NCL_PER_SUPERCL + tidxj+i;
264         ai = ci * CL_SIZE + tidxi;
265
266         xqbuf    = xq[ai] + (float4)(shift_vec[3 * nb_sci.shift], shift_vec[3 * nb_sci.shift + 1], shift_vec[3 * nb_sci.shift + 2], 0.0f);
267         xqbuf.w *= nbparam->epsfac;
268         xqib[(tidxj + i) * CL_SIZE + tidxi] = xqbuf;
269 #ifdef IATYPE_SHMEM
270 #ifndef LJ_COMB
271         /* Pre-load the i-atom types into shared memory */
272         atib[(tidxj + i) * CL_SIZE + tidxi]   = atom_types[ai];
273 #else
274         ljcpib[(tidxj + i) * CL_SIZE + tidxi] = lj_comb[ai];
275 #endif
276 #endif
277     }
278 #if !USE_SUBGROUP_ANY
279     /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
280     if (tidx == 0 || tidx == WARP_SIZE)
281     {
282         warp_any[widx] = 0;
283     }
284 #endif
285     barrier(CLK_LOCAL_MEM_FENCE);
286
287     for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
288     {
289         fci_buf[ci_offset] = (float3)(0.0f);
290     }
291
292 #ifdef LJ_EWALD
293     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
294     lje_coeff2   = nbparam->ewaldcoeff_lj*nbparam->ewaldcoeff_lj;
295     lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
296 #endif /* LJ_EWALD */
297
298
299 #ifdef CALC_ENERGIES
300     E_lj = 0.0f;
301     E_el = 0.0f;
302
303 #if defined EXCLUSION_FORCES /* Ewald or RF */
304     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
305     {
306         /* we have the diagonal: add the charge and LJ self interaction energy term */
307         for (i = 0; i < NCL_PER_SUPERCL; i++)
308         {
309 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
310             qi    = xqib[i * CL_SIZE + tidxi].w;
311             E_el += qi*qi;
312 #endif
313 #if defined LJ_EWALD
314             E_lj += nbfp_climg2d[atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2];
315 #endif      /* LJ_EWALD */
316         }
317
318         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
319 #ifdef LJ_EWALD
320         E_lj /= CL_SIZE;
321         E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
322 #endif  /* LJ_EWALD */
323
324 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
325         /* Correct for epsfac^2 due to adding qi^2 */
326         E_el /= nbparam->epsfac*CL_SIZE;
327 #if defined EL_RF || defined EL_CUTOFF
328         E_el *= -0.5f*c_rf;
329 #else
330         E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
331 #endif
332 #endif                                  /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
333     }
334 #endif                                  /* EXCLUSION_FORCES */
335
336 #endif                                  /* CALC_ENERGIES */
337
338 #ifdef EXCLUSION_FORCES
339     const int nonSelfInteraction  = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
340 #endif
341
342     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
343     for (j4 = cij4_start; j4 < cij4_end; j4++)
344     {
345         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
346         imask       = pl_cj4[j4].imei[widx].imask;
347         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
348
349         preloadCj4(&cjs, pl_cj4[j4].cj, tidxi, tidxj, imask != 0u);
350
351 #ifndef PRUNE_NBL
352         if (imask)
353 #endif
354         {
355             /* Unrolling this loop improves performance without pruning but
356              * with pruning it leads to slowdown.
357              *
358              * Tested with driver 1800.5
359              *
360              * TODO: check loop unrolling with NVIDIA OpenCL
361              */
362 #if !defined PRUNE_NBL && !defined _NVIDIA_SOURCE_
363 #pragma unroll 4
364 #endif
365             for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
366             {
367                 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
368                 {
369                     mask_ji = (1U << (jm * NCL_PER_SUPERCL));
370
371                     cj      = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
372                     aj      = cj * CL_SIZE + tidxj;
373
374                     /* load j atom data */
375                     xqbuf   = xq[aj];
376                     xj      = (float3)(xqbuf.xyz);
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 = (float3)(0.0f);
385
386 #if !defined PRUNE_NBL
387 #pragma unroll 8
388 #endif
389                     for (i = 0; i < NCL_PER_SUPERCL; i++)
390                     {
391                         if (imask & mask_ji)
392                         {
393                             ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
394                             ai      = ci * CL_SIZE + tidxi;      /* i atom index */
395
396                             /* all threads load an atom from i cluster ci into shmem! */
397                             xqbuf   = xqib[i * CL_SIZE + tidxi];
398                             xi      = (float3)(xqbuf.xyz);
399
400                             /* distance between i and j atoms */
401                             rv      = xi - xj;
402                             r2      = norm2(rv);
403
404 #ifdef PRUNE_NBL
405                             if (!gmx_sub_group_any(warp_any, widx, r2 < rlist_sq))
406                             {
407                                 imask &= ~mask_ji;
408                             }
409 #endif
410
411                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
412
413                             /* cutoff & exclusion check */
414 #ifdef EXCLUSION_FORCES
415                             if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
416 #else
417                             if ((r2 < rcoulomb_sq) * int_bit)
418 #endif
419                             {
420                                 /* load the rest of the i-atom parameters */
421                                 qi      = xqbuf.w;
422 #ifdef IATYPE_SHMEM
423 #ifndef LJ_COMB
424                                 typei   = atib[i * CL_SIZE + tidxi];
425 #else
426                                 ljcp_i  = ljcpib[i * CL_SIZE + tidxi];
427 #endif
428 #else                           /* IATYPE_SHMEM */
429 #ifndef LJ_COMB
430                                 typei   = atom_types[ai];
431 #else
432                                 ljcp_i  = lj_comb[ai];
433 #endif
434 #endif                          /* IATYPE_SHMEM */
435                                 /* LJ 6*C6 and 12*C12 */
436 #ifndef LJ_COMB
437                                 c6      = nbfp_climg2d[2 * (ntypes * typei + typej)];
438                                 c12     = nbfp_climg2d[2 * (ntypes * typei + typej)+1];
439 #else                           /* LJ_COMB */
440 #ifdef LJ_COMB_GEOM
441                                 c6      = ljcp_i.x * ljcp_j.x;
442                                 c12     = ljcp_i.y * ljcp_j.y;
443 #else
444                                 /* LJ 2^(1/6)*sigma and 12*epsilon */
445                                 sigma   = ljcp_i.x + ljcp_j.x;
446                                 epsilon = ljcp_i.y * ljcp_j.y;
447 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
448                                 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
449 #endif
450 #endif                          /* LJ_COMB_GEOM */
451 #endif                          /* LJ_COMB */
452
453                                 // Ensure distance do not become so small that r^-12 overflows
454                                 r2      = max(r2, NBNXN_MIN_RSQ);
455
456                                 inv_r   = rsqrt(r2);
457                                 inv_r2  = inv_r * inv_r;
458 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
459                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
460 #if defined EXCLUSION_FORCES
461                                 /* We could mask inv_r2, but with Ewald
462                                  * masking both inv_r6 and F_invr is faster */
463                                 inv_r6  *= int_bit;
464 #endif                          /* EXCLUSION_FORCES */
465
466                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
467 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
468                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot)*ONE_TWELVETH_F -
469                                                      c6 * (inv_r6 + nbparam->dispersion_shift.cpot)*ONE_SIXTH_F);
470
471 #endif
472 #else                           /* ! LJ_COMB_LB || CALC_ENERGIES */
473                                 float sig_r  = sigma*inv_r;
474                                 float sig_r2 = sig_r*sig_r;
475                                 float sig_r6 = sig_r2*sig_r2*sig_r2;
476 #if defined EXCLUSION_FORCES
477                                 sig_r6 *= int_bit;
478 #endif                          /* EXCLUSION_FORCES */
479
480                                 F_invr  = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
481 #endif                          /* ! LJ_COMB_LB || CALC_ENERGIES */
482
483
484 #ifdef LJ_FORCE_SWITCH
485 #ifdef CALC_ENERGIES
486                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
487 #else
488                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
489 #endif /* CALC_ENERGIES */
490 #endif /* LJ_FORCE_SWITCH */
491
492
493 #ifdef LJ_EWALD
494 #ifdef LJ_EWALD_COMB_GEOM
495 #ifdef CALC_ENERGIES
496                                 calculate_lj_ewald_comb_geom_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
497 #else
498                                 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
499 #endif                          /* CALC_ENERGIES */
500 #elif defined LJ_EWALD_COMB_LB
501                                 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
502 #ifdef CALC_ENERGIES
503                                                                int_bit, true, &F_invr, &E_lj_p
504 #else
505                                                                0, false, &F_invr, 0
506 #endif /* CALC_ENERGIES */
507                                                                );
508 #endif /* LJ_EWALD_COMB_GEOM */
509 #endif /* LJ_EWALD */
510
511 #ifdef LJ_POT_SWITCH
512 #ifdef CALC_ENERGIES
513                                 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
514 #else
515                                 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
516 #endif /* CALC_ENERGIES */
517 #endif /* LJ_POT_SWITCH */
518
519 #ifdef VDW_CUTOFF_CHECK
520                                 /* Separate VDW cut-off check to enable twin-range cut-offs
521                                  * (rvdw < rcoulomb <= rlist)
522                                  */
523                                 vdw_in_range  = (r2 < rvdw_sq) ? 1.0f : 0.0f;
524                                 F_invr       *= vdw_in_range;
525 #ifdef CALC_ENERGIES
526                                 E_lj_p       *= vdw_in_range;
527 #endif
528 #endif                          /* VDW_CUTOFF_CHECK */
529
530 #ifdef CALC_ENERGIES
531                                 E_lj    += E_lj_p;
532
533 #endif
534
535
536 #ifdef EL_CUTOFF
537 #ifdef EXCLUSION_FORCES
538                                 F_invr  += qi * qj_f * int_bit * inv_r2 * inv_r;
539 #else
540                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
541 #endif
542 #endif
543 #ifdef EL_RF
544                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
545 #endif
546 #if defined EL_EWALD_ANA
547                                 F_invr  += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
548 #elif defined EL_EWALD_TAB
549                                 F_invr  += qi * qj_f * (int_bit*inv_r2 -
550                                                         interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r, coulomb_tab_scale)
551                                                         ) * inv_r;
552 #endif                          /* EL_EWALD_ANA/TAB */
553
554 #ifdef CALC_ENERGIES
555 #ifdef EL_CUTOFF
556                                 E_el    += qi * qj_f * (int_bit*inv_r - c_rf);
557 #endif
558 #ifdef EL_RF
559                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
560 #endif
561 #ifdef EL_EWALD_ANY
562                                 /* 1.0f - erff is faster than erfcf */
563                                 E_el    += qi * qj_f * (inv_r * (int_bit - erf(r2 * inv_r * beta)) - int_bit * ewald_shift);
564 #endif                          /* EL_EWALD_ANY */
565 #endif
566                                 f_ij    = rv * F_invr;
567
568                                 /* accumulate j forces in registers */
569                                 fcj_buf -= f_ij;
570
571                                 /* accumulate i forces in registers */
572                                 fci_buf[i] += f_ij;
573                             }
574                         }
575
576                         /* shift the mask bit by 1 */
577                         mask_ji += mask_ji;
578                     }
579
580                     /* reduce j forces */
581                     reduce_force_j(f_buf, fcj_buf, f, tidxi, tidxj, aj);
582                 }
583             }
584 #ifdef PRUNE_NBL
585             /* Update the imask with the new one which does not contain the
586                out of range clusters anymore. */
587
588             pl_cj4[j4].imei[widx].imask = imask;
589 #endif
590         }
591     }
592
593     /* skip central shifts when summing shift forces */
594     if (nb_sci.shift == CENTRAL)
595     {
596         bCalcFshift = 0;
597     }
598     /* reduce i forces */
599     reduce_force_i_and_shift(f_buf, fci_buf, f, bCalcFshift != 0, tidxi, tidxj, sci,
600                              nb_sci.shift, fshift);
601
602 #ifdef CALC_ENERGIES
603     reduce_energy(f_buf, E_lj, E_el, e_lj, e_el, tidx);
604 #endif
605 }
606
607 #undef EL_EWALD_ANY
608 #undef EXCLUSION_FORCES
609 #undef LJ_EWALD
610
611 #undef LJ_COMB
612
613 #undef USE_CJ_PREFETCH