7df9231d515902ecb1f59b2df9302a15b4d6fe81
[alexxy/gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_cuda_kernel.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *  \brief
39  *  CUDA non-bonded kernel used through preprocessor-based code generation
40  *  of multiple kernel flavors, see nbnxn_cuda_kernels.cuh.
41  *
42  *  NOTE: No include fence as it is meant to be included multiple times.
43  *
44  *  \author Szilárd Páll <pall.szilard@gmail.com>
45  *  \author Berk Hess <hess@kth.se>
46  *  \ingroup module_nbnxm
47  */
48
49 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
50 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
51 #include "gromacs/gpu_utils/typecasts.cuh"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/pbcutil/ishift.h"
55 /* Note that floating-point constants in CUDA code should be suffixed
56  * with f (e.g. 0.5f), to stop the compiler producing intermediate
57  * code that is in double precision.
58  */
59
60 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
61 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
62 #    define EL_EWALD_ANY
63 #endif
64
65 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
66 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
67 #    define LJ_EWALD
68 #endif
69
70 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD \
71         || (defined EL_CUTOFF && defined CALC_ENERGIES)
72 /* Macro to control the calculation of exclusion forces in the kernel
73  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
74  * energy terms.
75  *
76  * Note: convenience macro, needs to be undef-ed at the end of the file.
77  */
78 #    define EXCLUSION_FORCES
79 #endif
80
81 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB
82 #    define LJ_COMB
83 #endif
84
85 /*
86    Kernel launch parameters:
87     - #blocks   = #pair lists, blockId = pair list Id
88     - #threads  = NTHREAD_Z * c_clSize^2
89     - shmem     = see nbnxn_cuda.cu:calc_shmem_required_nonbonded()
90
91     Each thread calculates an i force-component taking one pair of i-j atoms.
92  */
93
94 /**@{*/
95 /*! \brief Compute capability dependent definition of kernel launch configuration parameters.
96  *
97  * NTHREAD_Z controls the number of j-clusters processed concurrently on NTHREAD_Z
98  * warp-pairs per block.
99  *
100  * - On CC 3.0-3.5, and >=5.0 NTHREAD_Z == 1, translating to 64 th/block with 16
101  * blocks/multiproc, is the fastest even though this setup gives low occupancy
102  * (except on 6.0).
103  * NTHREAD_Z > 1 results in excessive register spilling unless the minimum blocks
104  * per multiprocessor is reduced proportionally to get the original number of max
105  * threads in flight (and slightly lower performance).
106  * - On CC 3.7 there are enough registers to double the number of threads; using
107  * NTHREADS_Z == 2 is fastest with 16 blocks (TODO: test with RF and other kernels
108  * with low-register use).
109  *
110  * Note that the current kernel implementation only supports NTHREAD_Z > 1 with
111  * shuffle-based reduction, hence CC >= 3.0.
112  *
113  *
114  * NOTEs on Volta / CUDA 9 extensions:
115  *
116  * - While active thread masks are required for the warp collectives
117  *   (we use any and shfl), the kernel is designed such that all conditions
118  *   (other than the inner-most distance check) including loop trip counts
119  *   are warp-synchronous. Therefore, we don't need ballot to compute the
120  *   active masks as these are all full-warp masks.
121  *
122  * - TODO: reconsider the use of __syncwarp(): its only role is currently to prevent
123  *   WAR hazard due to the cj preload; we should try to replace it with direct
124  *   loads (which may be faster given the improved L1 on Volta).
125  */
126
127 /* Kernel launch bounds for different compute capabilities. The value of NTHREAD_Z
128  * determines the number of threads per block and it is chosen such that
129  * 16 blocks/multiprocessor can be kept in flight.
130  * - CC 3.0,3.5, and >=5.0: NTHREAD_Z=1, (64, 16) bounds
131  * - CC 3.7:                NTHREAD_Z=2, (128, 16) bounds
132  *
133  * Note: convenience macros, need to be undef-ed at the end of the file.
134  */
135 #if GMX_PTX_ARCH == 370
136 #    define NTHREAD_Z (2)
137 #    define MIN_BLOCKS_PER_MP (16)
138 #else
139 #    define NTHREAD_Z (1)
140 #    define MIN_BLOCKS_PER_MP (16)
141 #endif /* GMX_PTX_ARCH == 370 */
142 #define THREADS_PER_BLOCK (c_clSize * c_clSize * NTHREAD_Z)
143
144 #if GMX_PTX_ARCH >= 350
145 /**@}*/
146 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
147 #else
148 __launch_bounds__(THREADS_PER_BLOCK)
149 #endif /* GMX_PTX_ARCH >= 350 */
150 #ifdef PRUNE_NBL
151 #    ifdef CALC_ENERGIES
152         __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
153 #    else
154         __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
155 #    endif /* CALC_ENERGIES */
156 #else
157 #    ifdef CALC_ENERGIES
158         __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
159 #    else
160         __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
161 #    endif /* CALC_ENERGIES */
162 #endif     /* PRUNE_NBL */
163                 (const NBAtomDataGpu atdat, const NBParamGpu nbparam, const Nbnxm::gpu_plist plist, bool bCalcFshift)
164 #ifdef FUNCTION_DECLARATION_ONLY
165                         ; /* Only do function declaration, omit the function body. */
166 #else
167 {
168     /* convenience variables */
169     const nbnxn_sci_t* pl_sci = plist.sci;
170 #    ifndef PRUNE_NBL
171     const
172 #    endif
173             nbnxn_cj4_t* pl_cj4      = plist.cj4;
174     const nbnxn_excl_t*  excl        = plist.excl;
175 #    ifndef LJ_COMB
176     const int*           atom_types  = atdat.atomTypes;
177     int                  ntypes      = atdat.numTypes;
178 #    else
179     const float2* lj_comb = atdat.ljComb;
180     float2        ljcp_i, ljcp_j;
181 #    endif
182     const float4*        xq          = atdat.xq;
183     float3*              f           = asFloat3(atdat.f);
184     const float3*        shift_vec   = asFloat3(atdat.shiftVec);
185     float                rcoulomb_sq = nbparam.rcoulomb_sq;
186 #    ifdef VDW_CUTOFF_CHECK
187     float                rvdw_sq     = nbparam.rvdw_sq;
188     float                vdw_in_range;
189 #    endif
190 #    ifdef LJ_EWALD
191     float                lje_coeff2, lje_coeff6_6;
192 #    endif
193 #    ifdef EL_RF
194     float                two_k_rf    = nbparam.two_k_rf;
195 #    endif
196 #    ifdef EL_EWALD_ANA
197     float                beta2       = nbparam.ewald_beta * nbparam.ewald_beta;
198     float                beta3       = nbparam.ewald_beta * nbparam.ewald_beta * nbparam.ewald_beta;
199 #    endif
200 #    ifdef PRUNE_NBL
201     float                rlist_sq    = nbparam.rlistOuter_sq;
202 #    endif
203
204 #    ifdef CALC_ENERGIES
205 #        ifdef EL_EWALD_ANY
206     float                beta        = nbparam.ewald_beta;
207     float                ewald_shift = nbparam.sh_ewald;
208 #        else
209     float reactionFieldShift = nbparam.c_rf;
210 #        endif /* EL_EWALD_ANY */
211     float*               e_lj        = atdat.eLJ;
212     float*               e_el        = atdat.eElec;
213 #    endif     /* CALC_ENERGIES */
214
215     /* thread/block/warp id-s */
216     unsigned int tidxi = threadIdx.x;
217     unsigned int tidxj = threadIdx.y;
218     unsigned int tidx  = threadIdx.y * blockDim.x + threadIdx.x;
219 #    if NTHREAD_Z == 1
220     unsigned int tidxz = 0;
221 #    else
222     unsigned int  tidxz = threadIdx.z;
223 #    endif
224     unsigned int bidx  = blockIdx.x;
225     unsigned int widx  = tidx / warp_size; /* warp index */
226
227     int          sci, ci, cj, ai, aj, cij4_start, cij4_end;
228 #    ifndef LJ_COMB
229     int          typei, typej;
230 #    endif
231     int          i, jm, j4, wexcl_idx;
232     float        qi, qj_f, r2, inv_r, inv_r2;
233 #    if !defined LJ_COMB_LB || defined CALC_ENERGIES
234     float        inv_r6, c6, c12;
235 #    endif
236 #    ifdef LJ_COMB_LB
237     float        sigma, epsilon;
238 #    endif
239     float        int_bit, F_invr;
240 #    ifdef CALC_ENERGIES
241     float        E_lj, E_el;
242 #    endif
243 #    if defined CALC_ENERGIES || defined LJ_POT_SWITCH
244     float        E_lj_p;
245 #    endif
246     unsigned int wexcl, imask, mask_ji;
247     float4       xqbuf;
248     float3       xi, xj, rv, f_ij, fcj_buf;
249     float3       fci_buf[c_nbnxnGpuNumClusterPerSupercluster]; /* i force buffer */
250     nbnxn_sci_t  nb_sci;
251
252     /*! i-cluster interaction mask for a super-cluster with all c_nbnxnGpuNumClusterPerSupercluster=8 bits set */
253     const unsigned superClInteractionMask = ((1U << c_nbnxnGpuNumClusterPerSupercluster) - 1U);
254
255     /*********************************************************************
256      * Set up shared memory pointers.
257      * sm_nextSlotPtr should always be updated to point to the "next slot",
258      * that is past the last point where data has been stored.
259      */
260     extern __shared__ char sm_dynamicShmem[];
261     char*                  sm_nextSlotPtr = sm_dynamicShmem;
262     static_assert(sizeof(char) == 1,
263                   "The shared memory offset calculation assumes that char is 1 byte");
264
265     /* shmem buffer for i x+q pre-loading */
266     float4* xqib = (float4*)sm_nextSlotPtr;
267     sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*xqib));
268
269     /* shmem buffer for cj, for each warp separately */
270     int* cjs = (int*)(sm_nextSlotPtr);
271     /* the cjs buffer's use expects a base pointer offset for pairs of warps in the j-concurrent execution */
272     cjs += tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
273     sm_nextSlotPtr += (NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(*cjs));
274
275 #    ifndef LJ_COMB
276     /* shmem buffer for i atom-type pre-loading */
277     int* atib = (int*)sm_nextSlotPtr;
278     sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*atib));
279 #    else
280     /* shmem buffer for i-atom LJ combination rule parameters */
281     float2* ljcpib = (float2*)sm_nextSlotPtr;
282     sm_nextSlotPtr += (c_nbnxnGpuNumClusterPerSupercluster * c_clSize * sizeof(*ljcpib));
283 #    endif
284     /*********************************************************************/
285
286     nb_sci     = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
287     sci        = nb_sci.sci;           /* super-cluster */
288     cij4_start = nb_sci.cj4_ind_start; /* first ...*/
289     cij4_end   = nb_sci.cj4_ind_end;   /* and last index of j clusters */
290
291     if (tidxz == 0)
292     {
293         /* Pre-load i-atom x and q into shared memory */
294         ci = sci * c_nbnxnGpuNumClusterPerSupercluster + tidxj;
295         ai = ci * c_clSize + tidxi;
296
297         float* shiftptr = (float*)&shift_vec[nb_sci.shift];
298         xqbuf = xq[ai] + make_float4(LDG(shiftptr), LDG(shiftptr + 1), LDG(shiftptr + 2), 0.0f);
299         xqbuf.w *= nbparam.epsfac;
300         xqib[tidxj * c_clSize + tidxi] = xqbuf;
301
302 #    ifndef LJ_COMB
303         /* Pre-load the i-atom types into shared memory */
304         atib[tidxj * c_clSize + tidxi] = atom_types[ai];
305 #    else
306         /* Pre-load the LJ combination parameters into shared memory */
307         ljcpib[tidxj * c_clSize + tidxi] = lj_comb[ai];
308 #    endif
309     }
310     __syncthreads();
311
312     for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
313     {
314         fci_buf[i] = make_float3(0.0f);
315     }
316
317 #    ifdef LJ_EWALD
318     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
319     lje_coeff2   = nbparam.ewaldcoeff_lj * nbparam.ewaldcoeff_lj;
320     lje_coeff6_6 = lje_coeff2 * lje_coeff2 * lje_coeff2 * c_oneSixth;
321 #    endif
322
323
324 #    ifdef CALC_ENERGIES
325     E_lj         = 0.0f;
326     E_el         = 0.0f;
327
328 #        ifdef EXCLUSION_FORCES /* Ewald or RF */
329     if (nb_sci.shift == gmx::c_centralShiftIndex
330         && pl_cj4[cij4_start].cj[0] == sci * c_nbnxnGpuNumClusterPerSupercluster)
331     {
332         /* we have the diagonal: add the charge and LJ self interaction energy term */
333         for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
334         {
335 #            if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
336             qi = xqib[i * c_clSize + tidxi].w;
337             E_el += qi * qi;
338 #            endif
339
340 #            ifdef LJ_EWALD
341             // load only the first 4 bytes of the parameter pair (equivalent with nbfp[idx].x)
342             E_lj += LDG((float*)&nbparam.nbfp[atom_types[(sci * c_nbnxnGpuNumClusterPerSupercluster + i) * c_clSize + tidxi]
343                                               * (ntypes + 1)]);
344 #            endif
345         }
346
347         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
348 #            ifdef LJ_EWALD
349         E_lj /= c_clSize * NTHREAD_Z;
350         E_lj *= 0.5f * c_oneSixth * lje_coeff6_6;
351 #            endif
352
353 #            if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
354         /* Correct for epsfac^2 due to adding qi^2 */
355         E_el /= nbparam.epsfac * c_clSize * NTHREAD_Z;
356 #                if defined EL_RF || defined EL_CUTOFF
357         E_el *= -0.5f * reactionFieldShift;
358 #                else
359         E_el *= -beta * M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
360 #                endif
361 #            endif /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
362     }
363 #        endif     /* EXCLUSION_FORCES */
364
365 #    endif /* CALC_ENERGIES */
366
367 #    ifdef EXCLUSION_FORCES
368     const int nonSelfInteraction = !(nb_sci.shift == gmx::c_centralShiftIndex & tidxj <= tidxi);
369 #    endif
370
371     /* loop over the j clusters = seen by any of the atoms in the current super-cluster;
372      * The loop stride NTHREAD_Z ensures that consecutive warps-pairs are assigned
373      * consecutive j4's entries.
374      */
375     for (j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
376     {
377         wexcl_idx = pl_cj4[j4].imei[widx].excl_ind;
378         imask     = pl_cj4[j4].imei[widx].imask;
379         wexcl     = excl[wexcl_idx].pair[(tidx) & (warp_size - 1)];
380
381 #    ifndef PRUNE_NBL
382         if (imask)
383 #    endif
384         {
385             /* Pre-load cj into shared memory on both warps separately */
386             if ((tidxj == 0 | tidxj == 4) & (tidxi < c_nbnxnGpuJgroupSize))
387             {
388                 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize / c_splitClSize] = pl_cj4[j4].cj[tidxi];
389             }
390             __syncwarp(c_fullWarpMask);
391
392             /* Unrolling this loop
393                - with pruning leads to register spilling;
394                - on Kepler and later it is much slower;
395                Tested with up to nvcc 7.5 */
396             for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
397             {
398                 if (imask & (superClInteractionMask << (jm * c_nbnxnGpuNumClusterPerSupercluster)))
399                 {
400                     mask_ji = (1U << (jm * c_nbnxnGpuNumClusterPerSupercluster));
401
402                     cj = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize / c_splitClSize];
403                     aj = cj * c_clSize + tidxj;
404
405                     /* load j atom data */
406                     xqbuf = xq[aj];
407                     xj    = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
408                     qj_f  = xqbuf.w;
409 #    ifndef LJ_COMB
410                     typej = atom_types[aj];
411 #    else
412                     ljcp_j = lj_comb[aj];
413 #    endif
414
415                     fcj_buf = make_float3(0.0f);
416
417 #    if !defined PRUNE_NBL
418 #        pragma unroll 8
419 #    endif
420                     for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
421                     {
422                         if (imask & mask_ji)
423                         {
424                             ci = sci * c_nbnxnGpuNumClusterPerSupercluster + i; /* i cluster index */
425
426                             /* all threads load an atom from i cluster ci into shmem! */
427                             xqbuf = xqib[i * c_clSize + tidxi];
428                             xi    = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
429
430                             /* distance between i and j atoms */
431                             rv = xi - xj;
432                             r2 = norm2(rv);
433
434 #    ifdef PRUNE_NBL
435                             /* If _none_ of the atoms pairs are in cutoff range,
436                                the bit corresponding to the current
437                                cluster-pair in imask gets set to 0. */
438                             if (!__any_sync(c_fullWarpMask, r2 < rlist_sq))
439                             {
440                                 imask &= ~mask_ji;
441                             }
442 #    endif
443
444                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
445
446                             /* cutoff & exclusion check */
447 #    ifdef EXCLUSION_FORCES
448                             if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
449 #    else
450                             if ((r2 < rcoulomb_sq) * int_bit)
451 #    endif
452                             {
453                                 /* load the rest of the i-atom parameters */
454                                 qi = xqbuf.w;
455
456 #    ifndef LJ_COMB
457                                 /* LJ 6*C6 and 12*C12 */
458                                 typei = atib[i * c_clSize + tidxi];
459                                 fetch_nbfp_c6_c12(c6, c12, nbparam, ntypes * typei + typej);
460 #    else
461                                 ljcp_i       = ljcpib[i * c_clSize + tidxi];
462 #        ifdef LJ_COMB_GEOM
463                                 c6           = ljcp_i.x * ljcp_j.x;
464                                 c12          = ljcp_i.y * ljcp_j.y;
465 #        else
466                                 /* LJ 2^(1/6)*sigma and 12*epsilon */
467                                 sigma   = ljcp_i.x + ljcp_j.x;
468                                 epsilon = ljcp_i.y * ljcp_j.y;
469 #            if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
470                                 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
471 #            endif
472 #        endif /* LJ_COMB_GEOM */
473 #    endif     /* LJ_COMB */
474
475                                 // Ensure distance do not become so small that r^-12 overflows
476                                 r2 = max(r2, c_nbnxnMinDistanceSquared);
477
478                                 inv_r  = rsqrt(r2);
479                                 inv_r2 = inv_r * inv_r;
480 #    if !defined LJ_COMB_LB || defined CALC_ENERGIES
481                                 inv_r6 = inv_r2 * inv_r2 * inv_r2;
482 #        ifdef EXCLUSION_FORCES
483                                 /* We could mask inv_r2, but with Ewald
484                                  * masking both inv_r6 and F_invr is faster */
485                                 inv_r6 *= int_bit;
486 #        endif /* EXCLUSION_FORCES */
487
488                                 F_invr = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
489 #        if defined CALC_ENERGIES || defined LJ_POT_SWITCH
490                                 E_lj_p = int_bit
491                                          * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot) * c_oneTwelveth
492                                             - c6 * (inv_r6 + nbparam.dispersion_shift.cpot) * c_oneSixth);
493 #        endif
494 #    else /* !LJ_COMB_LB || CALC_ENERGIES */
495                                 float sig_r  = sigma * inv_r;
496                                 float sig_r2 = sig_r * sig_r;
497                                 float sig_r6 = sig_r2 * sig_r2 * sig_r2;
498 #        ifdef EXCLUSION_FORCES
499                                 sig_r6 *= int_bit;
500 #        endif /* EXCLUSION_FORCES */
501
502                                 F_invr = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
503 #    endif     /* !LJ_COMB_LB || CALC_ENERGIES */
504
505 #    ifdef LJ_FORCE_SWITCH
506 #        ifdef CALC_ENERGIES
507                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
508 #        else
509                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
510 #        endif /* CALC_ENERGIES */
511 #    endif     /* LJ_FORCE_SWITCH */
512
513
514 #    ifdef LJ_EWALD
515 #        ifdef LJ_EWALD_COMB_GEOM
516 #            ifdef CALC_ENERGIES
517                                 calculate_lj_ewald_comb_geom_F_E(
518                                         nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
519 #            else
520                                 calculate_lj_ewald_comb_geom_F(
521                                         nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
522 #            endif /* CALC_ENERGIES */
523 #        elif defined LJ_EWALD_COMB_LB
524                                 calculate_lj_ewald_comb_LB_F_E(nbparam,
525                                                                typei,
526                                                                typej,
527                                                                r2,
528                                                                inv_r2,
529                                                                lje_coeff2,
530                                                                lje_coeff6_6,
531 #            ifdef CALC_ENERGIES
532                                                                int_bit,
533                                                                &F_invr,
534                                                                &E_lj_p
535 #            else
536                                                                0,
537                                                                &F_invr,
538                                                                nullptr
539 #            endif /* CALC_ENERGIES */
540                                 );
541 #        endif     /* LJ_EWALD_COMB_GEOM */
542 #    endif         /* LJ_EWALD */
543
544 #    ifdef LJ_POT_SWITCH
545 #        ifdef CALC_ENERGIES
546                                 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
547 #        else
548                                 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
549 #        endif /* CALC_ENERGIES */
550 #    endif     /* LJ_POT_SWITCH */
551
552 #    ifdef VDW_CUTOFF_CHECK
553                                 /* Separate VDW cut-off check to enable twin-range cut-offs
554                                  * (rvdw < rcoulomb <= rlist)
555                                  */
556                                 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
557                                 F_invr *= vdw_in_range;
558 #        ifdef CALC_ENERGIES
559                                 E_lj_p *= vdw_in_range;
560 #        endif
561 #    endif /* VDW_CUTOFF_CHECK */
562
563 #    ifdef CALC_ENERGIES
564                                 E_lj += E_lj_p;
565 #    endif
566
567
568 #    ifdef EL_CUTOFF
569 #        ifdef EXCLUSION_FORCES
570                                 F_invr += qi * qj_f * int_bit * inv_r2 * inv_r;
571 #        else
572                                 F_invr += qi * qj_f * inv_r2 * inv_r;
573 #        endif
574 #    endif
575 #    ifdef EL_RF
576                                 F_invr += qi * qj_f * (int_bit * inv_r2 * inv_r - two_k_rf);
577 #    endif
578 #    if defined   EL_EWALD_ANA
579                                 F_invr += qi * qj_f
580                                           * (int_bit * inv_r2 * inv_r + pmecorrF(beta2 * r2) * beta3);
581 #    elif defined EL_EWALD_TAB
582                                 F_invr += qi * qj_f
583                                           * (int_bit * inv_r2
584                                              - interpolate_coulomb_force_r(nbparam, r2 * inv_r))
585                                           * inv_r;
586 #    endif /* EL_EWALD_ANA/TAB */
587
588 #    ifdef CALC_ENERGIES
589 #        ifdef EL_CUTOFF
590                                 E_el += qi * qj_f * (int_bit * inv_r - reactionFieldShift);
591 #        endif
592 #        ifdef EL_RF
593                                 E_el += qi * qj_f
594                                         * (int_bit * inv_r + 0.5f * two_k_rf * r2 - reactionFieldShift);
595 #        endif
596 #        ifdef EL_EWALD_ANY
597                                 /* 1.0f - erff is faster than erfcf */
598                                 E_el += qi * qj_f
599                                         * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
600 #        endif /* EL_EWALD_ANY */
601 #    endif
602                                 f_ij = rv * F_invr;
603
604                                 /* accumulate j forces in registers */
605                                 fcj_buf -= f_ij;
606
607                                 /* accumulate i forces in registers */
608                                 fci_buf[i] += f_ij;
609                             }
610                         }
611
612                         /* shift the mask bit by 1 */
613                         mask_ji += mask_ji;
614                     }
615
616                     /* reduce j forces */
617                     reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj, c_fullWarpMask);
618                 }
619             }
620 #    ifdef PRUNE_NBL
621             /* Update the imask with the new one which does not contain the
622                out of range clusters anymore. */
623             pl_cj4[j4].imei[widx].imask = imask;
624 #    endif
625         }
626         // avoid shared memory WAR hazards between loop iterations
627         __syncwarp(c_fullWarpMask);
628     }
629
630     /* skip central shifts when summing shift forces */
631     if (nb_sci.shift == gmx::c_centralShiftIndex)
632     {
633         bCalcFshift = false;
634     }
635
636     float fshift_buf = 0.0f;
637
638     /* reduce i forces */
639     for (i = 0; i < c_nbnxnGpuNumClusterPerSupercluster; i++)
640     {
641         ai = (sci * c_nbnxnGpuNumClusterPerSupercluster + i) * c_clSize + tidxi;
642         reduce_force_i_warp_shfl(fci_buf[i], f, &fshift_buf, bCalcFshift, tidxj, ai, c_fullWarpMask);
643     }
644
645     /* add up local shift forces into global mem, tidxj indexes x,y,z */
646     if (bCalcFshift && (tidxj & 3) < 3)
647     {
648         float3* fShift = asFloat3(atdat.fShift);
649         atomicAdd(&(fShift[nb_sci.shift].x) + (tidxj & 3), fshift_buf);
650     }
651
652 #    ifdef CALC_ENERGIES
653     /* reduce the energies over warps and store into global memory */
654     reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx, c_fullWarpMask);
655 #    endif
656 }
657 #endif /* FUNCTION_DECLARATION_ONLY */
658
659 #undef NTHREAD_Z
660 #undef MIN_BLOCKS_PER_MP
661 #undef THREADS_PER_BLOCK
662
663 #undef EL_EWALD_ANY
664 #undef EXCLUSION_FORCES
665 #undef LJ_EWALD
666
667 #undef LJ_COMB