Refactor cj preload in the nonbonded OpenCL kernels
[gromacs.git] / src / gromacs / mdlib / nbnxn_ocl / nbnxn_ocl_kernel_amd.clh
blob22ba616cf01ae650f3281af01118ad44b4e0fc39
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  */
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  */
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 #define USE_CJ_PREFETCH 1
51 #include "nbnxn_ocl_kernel_utils.clh"
53 /////////////////////////////////////////////////////////////////////////////////////////////////
55 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
56 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
57 #define EL_EWALD_ANY
58 #endif
60 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
61 /* Macro to control the calculation of exclusion forces in the kernel
62  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
63  * energy terms.
64  *
65  * Note: convenience macro, needs to be undef-ed at the end of the file.
66  */
67 #define EXCLUSION_FORCES
68 #endif
70 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
71 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
72 #define LJ_EWALD
73 #endif
75 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
76 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
77 #define LJ_COMB
78 #endif
81    Kernel launch parameters:
82     - #blocks   = #pair lists, blockId = pair list Id
83     - #threads  = CL_SIZE^2
84     - shmem     = CL_SIZE^2 * sizeof(float)
86     Each thread calculates an i force-component taking one pair of i-j atoms.
88   TODO: implement 128 threads/wavefront by porting over the NTHREAD_Z/j4 loop
89   "horizontal splitting" over threads.
90  */
92 /* NOTE:
93  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.
94  Thus if more strings need to be appended a new macro must be written or it must be directly appended here.
96 __attribute__((reqd_work_group_size(CL_SIZE, CL_SIZE, 1)))
97 #ifdef PRUNE_NBL
98     #ifdef CALC_ENERGIES
99         __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_opencl)
100     #else
101         __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_opencl)
102     #endif
103 #else
104     #ifdef CALC_ENERGIES
105         __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_opencl)
106     #else
107         __kernel void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_opencl)
108     #endif
109 #endif
111 #ifndef LJ_COMB
112  int ntypes,                                                               /* IN  */
113 #endif
114  cl_nbparam_params_t nbparam_params,                                       /* IN  */
115  const __global float4 *restrict xq,                                       /* IN  */
116  __global float *restrict f,                /* stores float3 values */     /* OUT */
117  __global float *restrict e_lj,                                            /* OUT */
118  __global float *restrict e_el,                                            /* OUT */
119 __global float *restrict fshift,            /* stores float3 values */     /* OUT */
120 #ifdef LJ_COMB
121  const __global float2 *restrict lj_comb,      /* stores float2 values */  /* IN  */
122 #else
123  const __global int *restrict atom_types,                                  /* IN  */
124 #endif
125  const __global float *restrict shift_vec,  /* stores float3 values */     /* IN  */
126  __constant float* nbfp_climg2d,                                           /* IN  */
127  __constant float* nbfp_comb_climg2d,                                      /* IN  */
128  __constant float* coulomb_tab_climg2d,                                    /* IN  */
129  const __global nbnxn_sci_t* pl_sci,                                       /* IN  */
130 #ifndef PRUNE_NBL
131     const
132 #endif
133  __global nbnxn_cj4_t* pl_cj4,                                             /* OUT / IN */
134  const __global nbnxn_excl_t* excl,                                        /* IN  */
135  int bCalcFshift,                                                          /* IN  */
136  __local  float4   *xqib                                                   /* Pointer to dyn alloc'ed shmem */
139     /* convenience variables */
140     cl_nbparam_params_t *nbparam = &nbparam_params;
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
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 */
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;
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;
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);
216 #define LOCAL_OFFSET xqib + NCL_PER_SUPERCL * CL_SIZE
217     __local int *cjs;
218 #if USE_CJ_PREFETCH
219     /* shmem buffer for cj, for both warps separately */
220     cjs = (__local int *)(LOCAL_OFFSET);
221     #undef LOCAL_OFFSET
222     #define LOCAL_OFFSET cjs + 2 * NBNXN_GPU_JGROUP_SIZE
223 #endif
225 #ifdef IATYPE_SHMEM
226 #ifndef LJ_COMB
227     /* shmem buffer for i atom-type pre-loading */
228     __local int *atib      = (__local int *)(LOCAL_OFFSET);
229     #undef LOCAL_OFFSET
230     #define LOCAL_OFFSET atib + NCL_PER_SUPERCL * CL_SIZE
231 #else
232     __local float2 *ljcpib = (__local float2 *)(LOCAL_OFFSET);
233     #undef LOCAL_OFFSET
234     #define LOCAL_OFFSET ljcpib + NCL_PER_SUPERCL * CL_SIZE
235 #endif
236 #endif
238 #ifndef REDUCE_SHUFFLE
239     /* shmem j force buffer */
240     __local float *f_buf   = (__local float *)(LOCAL_OFFSET);
241     #undef LOCAL_OFFSET
242     #define LOCAL_OFFSET f_buf + CL_SIZE * CL_SIZE * 3
243 #endif
244     /* Local buffer used to implement __any warp vote function from CUDA.
245        volatile is used to avoid compiler optimizations for AMD builds. */
246     volatile __local uint *warp_any = (__local uint*)(LOCAL_OFFSET);
247 #undef LOCAL_OFFSET
249     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
250     sci         = nb_sci.sci;           /* super-cluster */
251     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
252     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
254     /* Pre-load i-atom x and q into shared memory */
255     ci = sci * NCL_PER_SUPERCL + tidxj;
256     ai = ci * CL_SIZE + tidxi;
258     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);
259     xqbuf.w *= nbparam->epsfac;
260     xqib[tidxj * CL_SIZE + tidxi] = xqbuf;
262 #ifdef IATYPE_SHMEM
263 #ifndef LJ_COMB
264     /* Pre-load the i-atom types into shared memory */
265     atib[tidxj * CL_SIZE + tidxi]   = atom_types[ai];
266 #else
267     ljcpib[tidxj * CL_SIZE + tidxi] = lj_comb[ai];
268 #endif
269 #endif
270     /* Initialise warp vote. (8x8 block) 2 warps for nvidia */
271     if(tidx==0 || tidx==32)
272         warp_any[widx] = 0;
274     barrier(CLK_LOCAL_MEM_FENCE);
276     for (ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
277     {
278         fci_buf[ci_offset] = (float3)(0.0f);
279     }
281 #ifdef LJ_EWALD
282     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
283     lje_coeff2   = nbparam->ewaldcoeff_lj*nbparam->ewaldcoeff_lj;
284     lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*ONE_SIXTH_F;
285 #endif /* LJ_EWALD */
288 #ifdef CALC_ENERGIES
289     E_lj = 0.0f;
290     E_el = 0.0f;
292 #if defined EXCLUSION_FORCES /* Ewald or RF */
293     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
294     {
295         /* we have the diagonal: add the charge and LJ self interaction energy term */
296         for (i = 0; i < NCL_PER_SUPERCL; i++)
297         {
298 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
299             qi    = xqib[i * CL_SIZE + tidxi].w;
300             E_el += qi*qi;
301 #endif
302 #if defined LJ_EWALD
303             E_lj += nbfp_climg2d[atom_types[(sci*NCL_PER_SUPERCL + i)*CL_SIZE + tidxi]*(ntypes + 1)*2];
304 #endif /* LJ_EWALD */
305         }
307         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
308 #ifdef LJ_EWALD
309         E_lj /= CL_SIZE;
310         E_lj *= 0.5f*ONE_SIXTH_F*lje_coeff6_6;
311 #endif  /* LJ_EWALD */
313 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
314         /* Correct for epsfac^2 due to adding qi^2 */
315         E_el /= nbparam->epsfac*CL_SIZE;
316 #if defined EL_RF || defined EL_CUTOFF
317         E_el *= -0.5f*c_rf;
318 #else
319         E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
320 #endif
321 #endif                                  /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
322     }
323 #endif                                  /* EXCLUSION_FORCES */
325 #endif                                  /* CALC_ENERGIES */
327 #ifdef EXCLUSION_FORCES
328     const int nonSelfInteraction  = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
329 #endif
331     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
332     for (j4 = cij4_start; j4 < cij4_end; j4++)
333     {
334         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
335         imask       = pl_cj4[j4].imei[widx].imask;
336         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
338         preloadCj4(cjs, pl_cj4[j4].cj, tidxi, tidxj, imask);
340 #ifndef PRUNE_NBL
341         if (imask)
342 #endif
343         {
344             /* Unrolling this loop improves performance without pruning but
345              * with pruning it leads to slowdown.
346              *
347              * Tested with driver 1800.5
348              */
349 #if !defined PRUNE_NBL
350 #pragma unroll 4
351 #endif
353             for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
354             {
355                 if (imask & (superClInteractionMask << (jm * NCL_PER_SUPERCL)))
356                 {
357                     mask_ji = (1U << (jm * NCL_PER_SUPERCL));
359                     cj      = loadCj(cjs, pl_cj4[j4].cj, jm, tidxi, tidxj);
360                     aj      = cj * CL_SIZE + tidxj;
362                     /* load j atom data */
363                     xqbuf   = xq[aj];
364                     xj      = (float3)(xqbuf.xyz);
365                     qj_f    = xqbuf.w;
366 #ifndef LJ_COMB
367                     typej   = atom_types[aj];
368 #else
369                     ljcp_j  = lj_comb[aj];
370 #endif
372                     fcj_buf = (float3)(0.0f);
374 #if !defined PRUNE_NBL
375 #pragma unroll 8
376 #endif
377                     for (i = 0; i < NCL_PER_SUPERCL; i++)
378                     {
379                         if (imask & mask_ji)
380                         {
381                             ci_offset   = i;                     /* i force buffer offset */
383                             ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
384                             ai      = ci * CL_SIZE + tidxi;      /* i atom index */
386                             /* all threads load an atom from i cluster ci into shmem! */
387                             xqbuf   = xqib[i * CL_SIZE + tidxi];
388                             xi      = (float3)(xqbuf.xyz);
390                             /* distance between i and j atoms */
391                             rv      = xi - xj;
392                             r2      = norm2(rv);
394 #ifdef PRUNE_NBL
395                             /* vote.. should code shmem serialisation, wonder what the hit will be */
396                             if (r2 < rlist_sq)
397                                 warp_any[widx]=1;
399                             /* If _none_ of the atoms pairs are in cutoff range,
400                                the bit corresponding to the current
401                                cluster-pair in imask gets set to 0. */
402                             if (!warp_any[widx])
403                                 imask &= ~mask_ji;
405                             warp_any[widx]=0;
406 #endif
408                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
410                             /* cutoff & exclusion check */
411 #ifdef EXCLUSION_FORCES
412                             if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
413 #else
414                             if ((r2 < rcoulomb_sq) * int_bit)
415 #endif
416                             {
417                                 /* load the rest of the i-atom parameters */
418                                 qi      = xqbuf.w;
419 #ifdef IATYPE_SHMEM
420 #ifndef LJ_COMB
421                                 typei   = atib[i * CL_SIZE + tidxi];
422 #else
423                                 ljcp_i  = ljcpib[i * CL_SIZE + tidxi];
424 #endif
425 #else /* IATYPE_SHMEM */
426 #ifndef LJ_COMB
427                                 typei   = atom_types[ai];
428 #else
429                                 ljcp_i  = lj_comb[ai];
430 #endif
431 #endif /* IATYPE_SHMEM */
432                                 /* LJ 6*C6 and 12*C12 */
433 #ifndef LJ_COMB
434                                 c6      = nbfp_climg2d[2 * (ntypes * typei + typej)];
435                                 c12     = nbfp_climg2d[2 * (ntypes * typei + typej)+1];
436 #else   /* LJ_COMB */
437 #ifdef LJ_COMB_GEOM
438                                 c6      = ljcp_i.x * ljcp_j.x;
439                                 c12     = ljcp_i.y * ljcp_j.y;
440 #else
441                                 /* LJ 2^(1/6)*sigma and 12*epsilon */
442                                 sigma   = ljcp_i.x + ljcp_j.x;
443                                 epsilon = ljcp_i.y * ljcp_j.y;
444 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
445                                 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
446 #endif
447 #endif                          /* LJ_COMB_GEOM */
448 #endif                          /* LJ_COMB */
450                                 // Ensure distance do not become so small that r^-12 overflows
451                                 r2      = max(r2,NBNXN_MIN_RSQ);
453                                 inv_r   = rsqrt(r2);
454                                 inv_r2  = inv_r * inv_r;
455 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
456                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
457 #if defined EXCLUSION_FORCES
458                                 /* We could mask inv_r2, but with Ewald
459                                  * masking both inv_r6 and F_invr is faster */
460                                 inv_r6  *= int_bit;
461 #endif                          /* EXCLUSION_FORCES */
463                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
464 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
465                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam->repulsion_shift.cpot)*ONE_TWELVETH_F -
466                                                      c6 * (inv_r6 + nbparam->dispersion_shift.cpot)*ONE_SIXTH_F);
468 #endif
469 #else                           /* ! LJ_COMB_LB || CALC_ENERGIES */
470                                 float sig_r  = sigma*inv_r;
471                                 float sig_r2 = sig_r*sig_r;
472                                 float sig_r6 = sig_r2*sig_r2*sig_r2;
473 #if defined EXCLUSION_FORCES
474                                 sig_r6 *= int_bit;
475 #endif                          /* EXCLUSION_FORCES */
477                                 F_invr  = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
478 #endif                          /* ! LJ_COMB_LB || CALC_ENERGIES */
481 #ifdef LJ_FORCE_SWITCH
482 #ifdef CALC_ENERGIES
483                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
484 #else
485                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
486 #endif /* CALC_ENERGIES */
487 #endif /* LJ_FORCE_SWITCH */
490 #ifdef LJ_EWALD
491 #ifdef LJ_EWALD_COMB_GEOM
492 #ifdef CALC_ENERGIES
493                                 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);
494 #else
495                                 calculate_lj_ewald_comb_geom_F(nbfp_comb_climg2d, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
496 #endif                          /* CALC_ENERGIES */
497 #elif defined LJ_EWALD_COMB_LB
498                                 calculate_lj_ewald_comb_LB_F_E(nbfp_comb_climg2d, nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
499 #ifdef CALC_ENERGIES
500                                                                int_bit, true, &F_invr, &E_lj_p
501 #else
502                                                                0, false, &F_invr, 0
503 #endif /* CALC_ENERGIES */
504                                                                );
505 #endif /* LJ_EWALD_COMB_GEOM */
506 #endif /* LJ_EWALD */
508 #ifdef LJ_POT_SWITCH
509 #ifdef CALC_ENERGIES
510                                 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
511 #else
512                                 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
513 #endif /* CALC_ENERGIES */
514 #endif /* LJ_POT_SWITCH */
516 #ifdef VDW_CUTOFF_CHECK
517                                 /* Separate VDW cut-off check to enable twin-range cut-offs
518                                  * (rvdw < rcoulomb <= rlist)
519                                  */
520                                 vdw_in_range  = (r2 < rvdw_sq) ? 1.0f : 0.0f;
521                                 F_invr       *= vdw_in_range;
522 #ifdef CALC_ENERGIES
523                                 E_lj_p       *= vdw_in_range;
524 #endif
525 #endif                          /* VDW_CUTOFF_CHECK */
527 #ifdef CALC_ENERGIES
528                                 E_lj    += E_lj_p;
530 #endif
533 #ifdef EL_CUTOFF
534 #ifdef EXCLUSION_FORCES
535                                 F_invr  += qi * qj_f * int_bit * inv_r2 * inv_r;
536 #else
537                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
538 #endif
539 #endif
540 #ifdef EL_RF
541                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
542 #endif
543 #if defined EL_EWALD_ANA
544                                 F_invr  += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
545 #elif defined EL_EWALD_TAB
546                                 F_invr  += qi * qj_f * (int_bit*inv_r2 -
547                                                         interpolate_coulomb_force_r(coulomb_tab_climg2d, r2 * inv_r, coulomb_tab_scale)
548                                                         ) * inv_r;
549 #endif /* EL_EWALD_ANA/TAB */
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;
565                                 /* accumulate j forces in registers */
566                                 fcj_buf -= f_ij;
568                                 /* accumulate i forces in registers */
569                                 fci_buf[ci_offset] += f_ij;
570                             }
571                         }
573                         /* shift the mask bit by 1 */
574                         mask_ji += mask_ji;
575                     }
577                     /* reduce j forces */
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;
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. */
591             pl_cj4[j4].imei[widx].imask = imask;
592 #endif
593         }
594     }
596     /* skip central shifts when summing shift forces */
597     if (nb_sci.shift == CENTRAL)
598     {
599         bCalcFshift = false;
600     }
602     fshift_buf = 0.0f;
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;
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     }
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     }
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);
636 #endif
639 #undef EL_EWALD_ANY
640 #undef EXCLUSION_FORCES
641 #undef LJ_EWALD
643 #undef LJ_COMB
645 #undef USE_CJ_PREFETCH