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