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