bf9cdf0d2eaed8bf5aa07736de68068005c17bc5
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_kernel.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,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
38  *  CUDA non-bonded kernel used through preprocessor-based code generation
39  *  of multiple kernel flavors, see nbnxn_cuda_kernels.cuh.
40  *
41  *  NOTE: No include fence as it is meant to be included multiple times.
42  *
43  *  \author Szilárd Páll <pall.szilard@gmail.com>
44  *  \author Berk Hess <hess@kth.se>
45  *  \ingroup module_mdlib
46  */
47
48 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
49 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/pbcutil/ishift.h"
52 /* Note that floating-point constants in CUDA code should be suffixed
53  * with f (e.g. 0.5f), to stop the compiler producing intermediate
54  * code that is in double precision.
55  */
56
57 #if GMX_PTX_ARCH < 300 && GMX_PTX_ARCH != 0
58 #error "nbnxn_cuda_kernel.cuh included with GMX_PTX_ARCH < 300 or host pass"
59 #endif
60
61 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
62 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
63 #define EL_EWALD_ANY
64 #endif
65
66 #if defined EL_EWALD_ANY || defined EL_RF || defined LJ_EWALD || (defined EL_CUTOFF && defined CALC_ENERGIES)
67 /* Macro to control the calculation of exclusion forces in the kernel
68  * We do that with Ewald (elec/vdw) and RF. Cut-off only has exclusion
69  * energy terms.
70  *
71  * Note: convenience macro, needs to be undef-ed at the end of the file.
72  */
73 #define EXCLUSION_FORCES
74 #endif
75
76 #if defined LJ_EWALD_COMB_GEOM || defined LJ_EWALD_COMB_LB
77 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
78 #define LJ_EWALD
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 2.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 #if (GMX_PTX_ARCH <= 210) && (NTHREAD_Z > 1)
146     #error NTHREAD_Z > 1 will give incorrect results on CC 2.x
147 #endif
148 /**@}*/
149 __launch_bounds__(THREADS_PER_BLOCK, MIN_BLOCKS_PER_MP)
150 #else
151 __launch_bounds__(THREADS_PER_BLOCK)
152 #endif /* GMX_PTX_ARCH >= 350 */
153 #ifdef PRUNE_NBL
154 #ifdef CALC_ENERGIES
155 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_prune_cuda)
156 #else
157 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_prune_cuda)
158 #endif /* CALC_ENERGIES */
159 #else
160 #ifdef CALC_ENERGIES
161 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _VF_cuda)
162 #else
163 __global__ void NB_KERNEL_FUNC_NAME(nbnxn_kernel, _F_cuda)
164 #endif /* CALC_ENERGIES */
165 #endif /* PRUNE_NBL */
166 (const cu_atomdata_t atdat,
167  const cu_nbparam_t nbparam,
168  const cu_plist_t plist,
169  bool bCalcFshift)
170 #ifdef FUNCTION_DECLARATION_ONLY
171 ;     /* Only do function declaration, omit the function body. */
172 #else
173 {
174     /* convenience variables */
175     const nbnxn_sci_t *pl_sci       = plist.sci;
176 #ifndef PRUNE_NBL
177     const
178 #endif
179     nbnxn_cj4_t        *pl_cj4      = plist.cj4;
180     const nbnxn_excl_t *excl        = plist.excl;
181 #ifndef LJ_COMB
182     const int          *atom_types  = atdat.atom_types;
183     int                 ntypes      = atdat.ntypes;
184 #else
185     const float2       *lj_comb     = atdat.lj_comb;
186     float2              ljcp_i, ljcp_j;
187 #endif
188     const float4       *xq          = atdat.xq;
189     float3             *f           = atdat.f;
190     const float3       *shift_vec   = atdat.shift_vec;
191     float               rcoulomb_sq = nbparam.rcoulomb_sq;
192 #ifdef VDW_CUTOFF_CHECK
193     float               rvdw_sq     = nbparam.rvdw_sq;
194     float               vdw_in_range;
195 #endif
196 #ifdef LJ_EWALD
197     float               lje_coeff2, lje_coeff6_6;
198 #endif
199 #ifdef EL_RF
200     float two_k_rf              = nbparam.two_k_rf;
201 #endif
202 #ifdef EL_EWALD_ANA
203     float beta2                 = nbparam.ewald_beta*nbparam.ewald_beta;
204     float beta3                 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
205 #endif
206 #ifdef PRUNE_NBL
207     float rlist_sq              = nbparam.rlistOuter_sq;
208 #endif
209
210 #ifdef CALC_ENERGIES
211 #ifdef EL_EWALD_ANY
212     float  beta        = nbparam.ewald_beta;
213     float  ewald_shift = nbparam.sh_ewald;
214 #else
215     float  c_rf        = nbparam.c_rf;
216 #endif /* EL_EWALD_ANY */
217     float *e_lj        = atdat.e_lj;
218     float *e_el        = atdat.e_el;
219 #endif /* CALC_ENERGIES */
220
221     /* thread/block/warp id-s */
222     unsigned int tidxi  = threadIdx.x;
223     unsigned int tidxj  = threadIdx.y;
224     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
225 #if NTHREAD_Z == 1
226     unsigned int tidxz  = 0;
227 #else
228     unsigned int tidxz  = threadIdx.z;
229 #endif
230     unsigned int bidx   = blockIdx.x;
231     unsigned int widx   = tidx / warp_size; /* warp index */
232
233     int          sci, ci, cj,
234                  ai, aj,
235                  cij4_start, cij4_end;
236 #ifndef LJ_COMB
237     int          typei, typej;
238 #endif
239     int          i, jm, j4, wexcl_idx;
240     float        qi, qj_f,
241                  r2, inv_r, inv_r2;
242 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
243     float        inv_r6, c6, c12;
244 #endif
245 #ifdef LJ_COMB_LB
246     float        sigma, epsilon;
247 #endif
248     float        int_bit,
249                  F_invr;
250 #ifdef CALC_ENERGIES
251     float        E_lj, E_el;
252 #endif
253 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
254     float        E_lj_p;
255 #endif
256     unsigned int wexcl, imask, mask_ji;
257     float4       xqbuf;
258     float3       xi, xj, rv, f_ij, fcj_buf;
259     float3       fci_buf[c_numClPerSupercl]; /* i force buffer */
260     nbnxn_sci_t  nb_sci;
261
262     /*! i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
263     const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
264
265     /*********************************************************************
266      * Set up shared memory pointers.
267      * sm_nextSlotPtr should always be updated to point to the "next slot",
268      * that is past the last point where data has been stored.
269      */
270     extern __shared__  char sm_dynamicShmem[];
271     char                   *sm_nextSlotPtr = sm_dynamicShmem;
272     static_assert(sizeof(char) == 1, "The shared memory offset calculation assumes that char is 1 byte");
273
274     /* shmem buffer for i x+q pre-loading */
275     float4 *xqib    = (float4 *)sm_nextSlotPtr;
276     sm_nextSlotPtr += (c_numClPerSupercl * c_clSize * sizeof(*xqib));
277
278     /* shmem buffer for cj, for each warp separately */
279     int *cjs        = (int *)(sm_nextSlotPtr);
280     /* the cjs buffer's use expects a base pointer offset for pairs of warps in the j-concurrent execution */
281     cjs            += tidxz * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize;
282     sm_nextSlotPtr += (NTHREAD_Z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(*cjs));
283
284 #ifndef LJ_COMB
285     /* shmem buffer for i atom-type pre-loading */
286     int *atib       = (int *)sm_nextSlotPtr;
287     sm_nextSlotPtr += (c_numClPerSupercl * c_clSize * sizeof(*atib));
288 #else
289     /* shmem buffer for i-atom LJ combination rule parameters */
290     float2 *ljcpib  = (float2 *)sm_nextSlotPtr;
291     sm_nextSlotPtr += (c_numClPerSupercl * c_clSize * sizeof(*ljcpib));
292 #endif
293     /*********************************************************************/
294
295     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
296     sci         = nb_sci.sci;           /* super-cluster */
297     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
298     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
299
300     if (tidxz == 0)
301     {
302         /* Pre-load i-atom x and q into shared memory */
303         ci = sci * c_numClPerSupercl + tidxj;
304         ai = ci * c_clSize + tidxi;
305
306         float  *shiftptr = (float *)&shift_vec[nb_sci.shift];
307         xqbuf    = xq[ai] + make_float4(LDG(shiftptr), LDG(shiftptr + 1), LDG(shiftptr + 2), 0.0f);
308         xqbuf.w *= nbparam.epsfac;
309         xqib[tidxj * c_clSize + tidxi] = xqbuf;
310
311 #ifndef LJ_COMB
312         /* Pre-load the i-atom types into shared memory */
313         atib[tidxj * c_clSize + tidxi] = atom_types[ai];
314 #else
315         /* Pre-load the LJ combination parameters into shared memory */
316         ljcpib[tidxj * c_clSize + tidxi] = lj_comb[ai];
317 #endif
318     }
319     __syncthreads();
320
321     for (i = 0; i < c_numClPerSupercl; i++)
322     {
323         fci_buf[i] = make_float3(0.0f);
324     }
325
326 #ifdef LJ_EWALD
327     /* TODO: we are trading registers with flops by keeping lje_coeff-s, try re-calculating it later */
328     lje_coeff2   = nbparam.ewaldcoeff_lj*nbparam.ewaldcoeff_lj;
329     lje_coeff6_6 = lje_coeff2*lje_coeff2*lje_coeff2*c_oneSixth;
330 #endif
331
332
333 #ifdef CALC_ENERGIES
334     E_lj = 0.0f;
335     E_el = 0.0f;
336
337 #ifdef EXCLUSION_FORCES /* Ewald or RF */
338     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*c_numClPerSupercl)
339     {
340         /* we have the diagonal: add the charge and LJ self interaction energy term */
341         for (i = 0; i < c_numClPerSupercl; i++)
342         {
343 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
344             qi    = xqib[i * c_clSize + tidxi].w;
345             E_el += qi*qi;
346 #endif
347
348 #ifdef LJ_EWALD
349     #if DISABLE_CUDA_TEXTURES
350             E_lj += LDG(&nbparam.nbfp[atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2]);
351     #else
352             E_lj += tex1Dfetch<float>(nbparam.nbfp_texobj, atom_types[(sci*c_numClPerSupercl + i)*c_clSize + tidxi]*(ntypes + 1)*2);
353     #endif
354 #endif
355         }
356
357         /* divide the self term(s) equally over the j-threads, then multiply with the coefficients. */
358 #ifdef LJ_EWALD
359         E_lj /= c_clSize*NTHREAD_Z;
360         E_lj *= 0.5f*c_oneSixth*lje_coeff6_6;
361 #endif
362
363 #if defined EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF
364         /* Correct for epsfac^2 due to adding qi^2 */
365         E_el /= nbparam.epsfac*c_clSize*NTHREAD_Z;
366 #if defined EL_RF || defined EL_CUTOFF
367         E_el *= -0.5f*c_rf;
368 #else
369         E_el *= -beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
370 #endif
371 #endif                                  /* EL_EWALD_ANY || defined EL_RF || defined EL_CUTOFF */
372     }
373 #endif                                  /* EXCLUSION_FORCES */
374
375 #endif                                  /* CALC_ENERGIES */
376
377 #ifdef EXCLUSION_FORCES
378     const int nonSelfInteraction = !(nb_sci.shift == CENTRAL & tidxj <= tidxi);
379 #endif
380
381     /* loop over the j clusters = seen by any of the atoms in the current super-cluster;
382      * The loop stride NTHREAD_Z ensures that consecutive warps-pairs are assigned
383      * consecutive j4's entries.
384      */
385     for (j4 = cij4_start + tidxz; j4 < cij4_end; j4 += NTHREAD_Z)
386     {
387         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
388         imask       = pl_cj4[j4].imei[widx].imask;
389         wexcl       = excl[wexcl_idx].pair[(tidx) & (warp_size - 1)];
390
391 #ifndef PRUNE_NBL
392         if (imask)
393 #endif
394         {
395             /* Pre-load cj into shared memory on both warps separately */
396             if ((tidxj == 0 | tidxj == 4) & (tidxi < c_nbnxnGpuJgroupSize))
397             {
398                 cjs[tidxi + tidxj * c_nbnxnGpuJgroupSize/c_splitClSize] = pl_cj4[j4].cj[tidxi];
399             }
400             gmx_syncwarp(c_fullWarpMask);
401
402             /* Unrolling this loop
403                - with pruning leads to register spilling;
404                - on Kepler and later it is much slower;
405                Tested with up to nvcc 7.5 */
406             for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
407             {
408                 if (imask & (superClInteractionMask << (jm * c_numClPerSupercl)))
409                 {
410                     mask_ji = (1U << (jm * c_numClPerSupercl));
411
412                     cj      = cjs[jm + (tidxj & 4) * c_nbnxnGpuJgroupSize/c_splitClSize];
413                     aj      = cj * c_clSize + tidxj;
414
415                     /* load j atom data */
416                     xqbuf   = xq[aj];
417                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
418                     qj_f    = xqbuf.w;
419 #ifndef LJ_COMB
420                     typej   = atom_types[aj];
421 #else
422                     ljcp_j  = lj_comb[aj];
423 #endif
424
425                     fcj_buf = make_float3(0.0f);
426
427 #if !defined PRUNE_NBL
428 #pragma unroll 8
429 #endif
430                     for (i = 0; i < c_numClPerSupercl; i++)
431                     {
432                         if (imask & mask_ji)
433                         {
434                             ci      = sci * c_numClPerSupercl + i; /* i cluster index */
435
436                             /* all threads load an atom from i cluster ci into shmem! */
437                             xqbuf   = xqib[i * c_clSize + tidxi];
438                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
439
440                             /* distance between i and j atoms */
441                             rv      = xi - xj;
442                             r2      = norm2(rv);
443
444 #ifdef PRUNE_NBL
445                             /* If _none_ of the atoms pairs are in cutoff range,
446                                the bit corresponding to the current
447                                cluster-pair in imask gets set to 0. */
448                             if (!gmx_any_sync(c_fullWarpMask, r2 < rlist_sq))
449                             {
450                                 imask &= ~mask_ji;
451                             }
452 #endif
453
454                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
455
456                             /* cutoff & exclusion check */
457 #ifdef EXCLUSION_FORCES
458                             if ((r2 < rcoulomb_sq) * (nonSelfInteraction | (ci != cj)))
459 #else
460                             if ((r2 < rcoulomb_sq) * int_bit)
461 #endif
462                             {
463                                 /* load the rest of the i-atom parameters */
464                                 qi      = xqbuf.w;
465
466 #ifndef LJ_COMB
467                                 /* LJ 6*C6 and 12*C12 */
468                                 typei   = atib[i * c_clSize + tidxi];
469                                 fetch_nbfp_c6_c12(c6, c12, nbparam, ntypes * typei + typej);
470 #else
471                                 ljcp_i  = ljcpib[i * c_clSize + tidxi];
472 #ifdef LJ_COMB_GEOM
473                                 c6      = ljcp_i.x * ljcp_j.x;
474                                 c12     = ljcp_i.y * ljcp_j.y;
475 #else
476                                 /* LJ 2^(1/6)*sigma and 12*epsilon */
477                                 sigma   = ljcp_i.x + ljcp_j.x;
478                                 epsilon = ljcp_i.y * ljcp_j.y;
479 #if defined CALC_ENERGIES || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
480                                 convert_sigma_epsilon_to_c6_c12(sigma, epsilon, &c6, &c12);
481 #endif
482 #endif                          /* LJ_COMB_GEOM */
483 #endif                          /* LJ_COMB */
484
485                                 // Ensure distance do not become so small that r^-12 overflows
486                                 r2      = max(r2, NBNXN_MIN_RSQ);
487
488                                 inv_r   = rsqrt(r2);
489                                 inv_r2  = inv_r * inv_r;
490 #if !defined LJ_COMB_LB || defined CALC_ENERGIES
491                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
492 #ifdef EXCLUSION_FORCES
493                                 /* We could mask inv_r2, but with Ewald
494                                  * masking both inv_r6 and F_invr is faster */
495                                 inv_r6  *= int_bit;
496 #endif                          /* EXCLUSION_FORCES */
497
498                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
499 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
500                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 + nbparam.repulsion_shift.cpot)*c_oneTwelveth -
501                                                      c6 * (inv_r6 + nbparam.dispersion_shift.cpot)*c_oneSixth);
502 #endif
503 #else                           /* !LJ_COMB_LB || CALC_ENERGIES */
504                                 float sig_r  = sigma*inv_r;
505                                 float sig_r2 = sig_r*sig_r;
506                                 float sig_r6 = sig_r2*sig_r2*sig_r2;
507 #ifdef EXCLUSION_FORCES
508                                 sig_r6 *= int_bit;
509 #endif                          /* EXCLUSION_FORCES */
510
511                                 F_invr  = epsilon * sig_r6 * (sig_r6 - 1.0f) * inv_r2;
512 #endif                          /* !LJ_COMB_LB || CALC_ENERGIES */
513
514 #ifdef LJ_FORCE_SWITCH
515 #ifdef CALC_ENERGIES
516                                 calculate_force_switch_F_E(nbparam, c6, c12, inv_r, r2, &F_invr, &E_lj_p);
517 #else
518                                 calculate_force_switch_F(nbparam, c6, c12, inv_r, r2, &F_invr);
519 #endif /* CALC_ENERGIES */
520 #endif /* LJ_FORCE_SWITCH */
521
522
523 #ifdef LJ_EWALD
524 #ifdef LJ_EWALD_COMB_GEOM
525 #ifdef CALC_ENERGIES
526                                 calculate_lj_ewald_comb_geom_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, int_bit, &F_invr, &E_lj_p);
527 #else
528                                 calculate_lj_ewald_comb_geom_F(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6, &F_invr);
529 #endif                          /* CALC_ENERGIES */
530 #elif defined LJ_EWALD_COMB_LB
531                                 calculate_lj_ewald_comb_LB_F_E(nbparam, typei, typej, r2, inv_r2, lje_coeff2, lje_coeff6_6,
532 #ifdef CALC_ENERGIES
533                                                                int_bit, &F_invr, &E_lj_p
534 #else
535                                                                0, &F_invr, nullptr
536 #endif /* CALC_ENERGIES */
537                                                                );
538 #endif /* LJ_EWALD_COMB_GEOM */
539 #endif /* LJ_EWALD */
540
541 #ifdef LJ_POT_SWITCH
542 #ifdef CALC_ENERGIES
543                                 calculate_potential_switch_F_E(nbparam, inv_r, r2, &F_invr, &E_lj_p);
544 #else
545                                 calculate_potential_switch_F(nbparam, inv_r, r2, &F_invr, &E_lj_p);
546 #endif /* CALC_ENERGIES */
547 #endif /* LJ_POT_SWITCH */
548
549 #ifdef VDW_CUTOFF_CHECK
550                                 /* Separate VDW cut-off check to enable twin-range cut-offs
551                                  * (rvdw < rcoulomb <= rlist)
552                                  */
553                                 vdw_in_range  = (r2 < rvdw_sq) ? 1.0f : 0.0f;
554                                 F_invr       *= vdw_in_range;
555 #ifdef CALC_ENERGIES
556                                 E_lj_p       *= vdw_in_range;
557 #endif
558 #endif                          /* VDW_CUTOFF_CHECK */
559
560 #ifdef CALC_ENERGIES
561                                 E_lj    += E_lj_p;
562 #endif
563
564
565 #ifdef EL_CUTOFF
566 #ifdef EXCLUSION_FORCES
567                                 F_invr  += qi * qj_f * int_bit * inv_r2 * inv_r;
568 #else
569                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
570 #endif
571 #endif
572 #ifdef EL_RF
573                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
574 #endif
575 #if defined EL_EWALD_ANA
576                                 F_invr  += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
577 #elif defined EL_EWALD_TAB
578                                 F_invr  += qi * qj_f * (int_bit*inv_r2 -
579                                                         interpolate_coulomb_force_r(nbparam, r2 * inv_r)) * inv_r;
580 #endif                          /* EL_EWALD_ANA/TAB */
581
582 #ifdef CALC_ENERGIES
583 #ifdef EL_CUTOFF
584                                 E_el    += qi * qj_f * (int_bit*inv_r - c_rf);
585 #endif
586 #ifdef EL_RF
587                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
588 #endif
589 #ifdef EL_EWALD_ANY
590                                 /* 1.0f - erff is faster than erfcf */
591                                 E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
592 #endif                          /* EL_EWALD_ANY */
593 #endif
594                                 f_ij    = rv * F_invr;
595
596                                 /* accumulate j forces in registers */
597                                 fcj_buf -= f_ij;
598
599                                 /* accumulate i forces in registers */
600                                 fci_buf[i] += f_ij;
601                             }
602                         }
603
604                         /* shift the mask bit by 1 */
605                         mask_ji += mask_ji;
606                     }
607
608                     /* reduce j forces */
609                     reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj, c_fullWarpMask);
610                 }
611             }
612 #ifdef PRUNE_NBL
613             /* Update the imask with the new one which does not contain the
614                out of range clusters anymore. */
615             pl_cj4[j4].imei[widx].imask = imask;
616 #endif
617         }
618         // avoid shared memory WAR hazards between loop iterations
619         gmx_syncwarp(c_fullWarpMask);
620     }
621
622     /* skip central shifts when summing shift forces */
623     if (nb_sci.shift == CENTRAL)
624     {
625         bCalcFshift = false;
626     }
627
628     float fshift_buf = 0.0f;
629
630     /* reduce i forces */
631     for (i = 0; i < c_numClPerSupercl; i++)
632     {
633         ai  = (sci * c_numClPerSupercl + i) * c_clSize + tidxi;
634         reduce_force_i_warp_shfl(fci_buf[i], f,
635                                  &fshift_buf, bCalcFshift,
636                                  tidxj, ai, c_fullWarpMask);
637     }
638
639     /* add up local shift forces into global mem, tidxj indexes x,y,z */
640     if (bCalcFshift && (tidxj & 3) < 3)
641     {
642         atomicAdd(&(atdat.fshift[nb_sci.shift].x) + (tidxj & 3), fshift_buf);
643     }
644
645 #ifdef CALC_ENERGIES
646     /* reduce the energies over warps and store into global memory */
647     reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx, c_fullWarpMask);
648 #endif
649 }
650 #endif /* FUNCTION_DECLARATION_ONLY */
651
652 #undef NTHREAD_Z
653 #undef MIN_BLOCKS_PER_MP
654 #undef THREADS_PER_BLOCK
655
656 #undef EL_EWALD_ANY
657 #undef EXCLUSION_FORCES
658 #undef LJ_EWALD
659
660 #undef LJ_COMB