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