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