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