2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 * Utility constant and function declaration for the CUDA non-bonded kernels.
39 * This header should be included once at the top level, just before the
40 * kernels are included (has to be preceded by nbnxn_cuda_types.h).
42 * \author Szilárd Páll <pall.szilard@gmail.com>
43 * \ingroup module_mdlib
47 /* Note that floating-point constants in CUDA code should be suffixed
48 * with f (e.g. 0.5f), to stop the compiler producing intermediate
49 * code that is in double precision.
52 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
53 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
54 #include "gromacs/gpu_utils/vectype_ops.cuh"
56 #include "nbnxn_cuda_types.h"
58 #ifndef NBNXN_CUDA_KERNEL_UTILS_CUH
59 #define NBNXN_CUDA_KERNEL_UTILS_CUH
61 /*! \brief Log of the i and j cluster size.
62 * change this together with c_clSize !*/
63 static const int c_clSizeLog2 = 3;
64 /*! \brief Square of cluster size. */
65 static const int c_clSizeSq = c_clSize*c_clSize;
66 /*! \brief j-cluster size after split (4 in the current implementation). */
67 static const int c_splitClSize = c_clSize/c_nbnxnGpuClusterpairSplit;
68 /*! \brief Stride in the force accumualation buffer */
69 static const int c_fbufStride = c_clSizeSq;
70 /*! \brief i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
71 static const unsigned superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
73 static const float c_oneSixth = 0.16666667f;
74 static const float c_oneTwelveth = 0.08333333f;
77 /*! Convert LJ sigma,epsilon parameters to C6,C12. */
78 static __forceinline__ __device__
79 void convert_sigma_epsilon_to_c6_c12(const float sigma,
86 sigma2 = sigma * sigma;
87 sigma6 = sigma2 *sigma2 * sigma2;
88 *c6 = epsilon * sigma6;
92 /*! Apply force switch, force + energy version. */
93 static __forceinline__ __device__
94 void calculate_force_switch_F(const cu_nbparam_t nbparam,
103 /* force switch constants */
104 float disp_shift_V2 = nbparam.dispersion_shift.c2;
105 float disp_shift_V3 = nbparam.dispersion_shift.c3;
106 float repu_shift_V2 = nbparam.repulsion_shift.c2;
107 float repu_shift_V3 = nbparam.repulsion_shift.c3;
110 r_switch = r - nbparam.rvdw_switch;
111 r_switch = r_switch >= 0.0f ? r_switch : 0.0f;
114 -c6*(disp_shift_V2 + disp_shift_V3*r_switch)*r_switch*r_switch*inv_r +
115 c12*(-repu_shift_V2 + repu_shift_V3*r_switch)*r_switch*r_switch*inv_r;
118 /*! Apply force switch, force-only version. */
119 static __forceinline__ __device__
120 void calculate_force_switch_F_E(const cu_nbparam_t nbparam,
130 /* force switch constants */
131 float disp_shift_V2 = nbparam.dispersion_shift.c2;
132 float disp_shift_V3 = nbparam.dispersion_shift.c3;
133 float repu_shift_V2 = nbparam.repulsion_shift.c2;
134 float repu_shift_V3 = nbparam.repulsion_shift.c3;
136 float disp_shift_F2 = nbparam.dispersion_shift.c2/3;
137 float disp_shift_F3 = nbparam.dispersion_shift.c3/4;
138 float repu_shift_F2 = nbparam.repulsion_shift.c2/3;
139 float repu_shift_F3 = nbparam.repulsion_shift.c3/4;
142 r_switch = r - nbparam.rvdw_switch;
143 r_switch = r_switch >= 0.0f ? r_switch : 0.0f;
146 -c6*(disp_shift_V2 + disp_shift_V3*r_switch)*r_switch*r_switch*inv_r +
147 c12*(-repu_shift_V2 + repu_shift_V3*r_switch)*r_switch*r_switch*inv_r;
149 c6*(disp_shift_F2 + disp_shift_F3*r_switch)*r_switch*r_switch*r_switch -
150 c12*(repu_shift_F2 + repu_shift_F3*r_switch)*r_switch*r_switch*r_switch;
153 /*! Apply potential switch, force-only version. */
154 static __forceinline__ __device__
155 void calculate_potential_switch_F(const cu_nbparam_t nbparam,
164 /* potential switch constants */
165 float switch_V3 = nbparam.vdw_switch.c3;
166 float switch_V4 = nbparam.vdw_switch.c4;
167 float switch_V5 = nbparam.vdw_switch.c5;
168 float switch_F2 = 3*nbparam.vdw_switch.c3;
169 float switch_F3 = 4*nbparam.vdw_switch.c4;
170 float switch_F4 = 5*nbparam.vdw_switch.c5;
173 r_switch = r - nbparam.rvdw_switch;
175 /* Unlike in the F+E kernel, conditional is faster here */
178 sw = 1.0f + (switch_V3 + (switch_V4 + switch_V5*r_switch)*r_switch)*r_switch*r_switch*r_switch;
179 dsw = (switch_F2 + (switch_F3 + switch_F4*r_switch)*r_switch)*r_switch*r_switch;
181 *F_invr = (*F_invr)*sw - inv_r*(*E_lj)*dsw;
185 /*! Apply potential switch, force + energy version. */
186 static __forceinline__ __device__
187 void calculate_potential_switch_F_E(const cu_nbparam_t nbparam,
196 /* potential switch constants */
197 float switch_V3 = nbparam.vdw_switch.c3;
198 float switch_V4 = nbparam.vdw_switch.c4;
199 float switch_V5 = nbparam.vdw_switch.c5;
200 float switch_F2 = 3*nbparam.vdw_switch.c3;
201 float switch_F3 = 4*nbparam.vdw_switch.c4;
202 float switch_F4 = 5*nbparam.vdw_switch.c5;
205 r_switch = r - nbparam.rvdw_switch;
206 r_switch = r_switch >= 0.0f ? r_switch : 0.0f;
208 /* Unlike in the F-only kernel, masking is faster here */
209 sw = 1.0f + (switch_V3 + (switch_V4 + switch_V5*r_switch)*r_switch)*r_switch*r_switch*r_switch;
210 dsw = (switch_F2 + (switch_F3 + switch_F4*r_switch)*r_switch)*r_switch*r_switch;
212 *F_invr = (*F_invr)*sw - inv_r*(*E_lj)*dsw;
217 /*! \brief Fetch C6 grid contribution coefficients and return the product of these.
219 * Depending on what is supported, it fetches parameters either
220 * using direct load, texture objects, or texrefs.
222 static __forceinline__ __device__
223 float calculate_lj_ewald_c6grid(const cu_nbparam_t nbparam,
227 #if DISABLE_CUDA_TEXTURES
228 return LDG(&nbparam.nbfp_comb[2*typei]) * LDG(&nbparam.nbfp_comb[2*typej]);
230 return tex1Dfetch<float>(nbparam.nbfp_comb_texobj, 2*typei) * tex1Dfetch<float>(nbparam.nbfp_comb_texobj, 2*typej);
231 #endif /* DISABLE_CUDA_TEXTURES */
235 /*! Calculate LJ-PME grid force contribution with
236 * geometric combination rule.
238 static __forceinline__ __device__
239 void calculate_lj_ewald_comb_geom_F(const cu_nbparam_t nbparam,
248 float c6grid, inv_r6_nm, cr2, expmcr2, poly;
250 c6grid = calculate_lj_ewald_c6grid(nbparam, typei, typej);
252 /* Recalculate inv_r6 without exclusion mask */
253 inv_r6_nm = inv_r2*inv_r2*inv_r2;
255 expmcr2 = expf(-cr2);
256 poly = 1.0f + cr2 + 0.5f*cr2*cr2;
258 /* Subtract the grid force from the total LJ force */
259 *F_invr += c6grid*(inv_r6_nm - expmcr2*(inv_r6_nm*poly + lje_coeff6_6))*inv_r2;
263 /*! Calculate LJ-PME grid force + energy contribution with
264 * geometric combination rule.
266 static __forceinline__ __device__
267 void calculate_lj_ewald_comb_geom_F_E(const cu_nbparam_t nbparam,
278 float c6grid, inv_r6_nm, cr2, expmcr2, poly, sh_mask;
280 c6grid = calculate_lj_ewald_c6grid(nbparam, typei, typej);
282 /* Recalculate inv_r6 without exclusion mask */
283 inv_r6_nm = inv_r2*inv_r2*inv_r2;
285 expmcr2 = expf(-cr2);
286 poly = 1.0f + cr2 + 0.5f*cr2*cr2;
288 /* Subtract the grid force from the total LJ force */
289 *F_invr += c6grid*(inv_r6_nm - expmcr2*(inv_r6_nm*poly + lje_coeff6_6))*inv_r2;
291 /* Shift should be applied only to real LJ pairs */
292 sh_mask = nbparam.sh_lj_ewald*int_bit;
293 *E_lj += c_oneSixth*c6grid*(inv_r6_nm*(1.0f - expmcr2*poly) + sh_mask);
296 /*! Fetch per-type LJ parameters.
298 * Depending on what is supported, it fetches parameters either
299 * using direct load, texture objects, or texrefs.
301 static __forceinline__ __device__
302 float2 fetch_nbfp_comb_c6_c12(const cu_nbparam_t nbparam,
306 #if DISABLE_CUDA_TEXTURES
307 /* Force an 8-byte fetch to save a memory instruction. */
308 float2 *nbfp_comb = (float2 *)nbparam.nbfp_comb;
309 c6c12 = LDG(&nbfp_comb[type]);
311 /* NOTE: as we always do 8-byte aligned loads, we could
312 fetch float2 here too just as above. */
313 c6c12.x = tex1Dfetch<float>(nbparam.nbfp_comb_texobj, 2*type);
314 c6c12.y = tex1Dfetch<float>(nbparam.nbfp_comb_texobj, 2*type + 1);
315 #endif /* DISABLE_CUDA_TEXTURES */
321 /*! Calculate LJ-PME grid force + energy contribution (if E_lj != NULL) with
322 * Lorentz-Berthelot combination rule.
323 * We use a single F+E kernel with conditional because the performance impact
324 * of this is pretty small and LB on the CPU is anyway very slow.
326 static __forceinline__ __device__
327 void calculate_lj_ewald_comb_LB_F_E(const cu_nbparam_t nbparam,
338 float c6grid, inv_r6_nm, cr2, expmcr2, poly;
339 float sigma, sigma2, epsilon;
341 /* sigma and epsilon are scaled to give 6*C6 */
342 float2 c6c12_i = fetch_nbfp_comb_c6_c12(nbparam, typei);
343 float2 c6c12_j = fetch_nbfp_comb_c6_c12(nbparam, typej);
345 sigma = c6c12_i.x + c6c12_j.x;
346 epsilon = c6c12_i.y * c6c12_j.y;
348 sigma2 = sigma*sigma;
349 c6grid = epsilon*sigma2*sigma2*sigma2;
351 /* Recalculate inv_r6 without exclusion mask */
352 inv_r6_nm = inv_r2*inv_r2*inv_r2;
354 expmcr2 = expf(-cr2);
355 poly = 1.0f + cr2 + 0.5f*cr2*cr2;
357 /* Subtract the grid force from the total LJ force */
358 *F_invr += c6grid*(inv_r6_nm - expmcr2*(inv_r6_nm*poly + lje_coeff6_6))*inv_r2;
364 /* Shift should be applied only to real LJ pairs */
365 sh_mask = nbparam.sh_lj_ewald*int_bit;
366 *E_lj += c_oneSixth*c6grid*(inv_r6_nm*(1.0f - expmcr2*poly) + sh_mask);
371 /*! Fetch two consecutive values from the Ewald correction F*r table.
373 * Depending on what is supported, it fetches parameters either
374 * using direct load, texture objects, or texrefs.
376 static __forceinline__ __device__
377 float2 fetch_coulomb_force_r(const cu_nbparam_t nbparam,
382 #if DISABLE_CUDA_TEXTURES
383 /* Can't do 8-byte fetch because some of the addresses will be misaligned. */
384 d.x = LDG(&nbparam.coulomb_tab[index]);
385 d.y = LDG(&nbparam.coulomb_tab[index + 1]);
387 d.x = tex1Dfetch<float>(nbparam.coulomb_tab_texobj, index);
388 d.y = tex1Dfetch<float>(nbparam.coulomb_tab_texobj, index + 1);
389 #endif // DISABLE_CUDA_TEXTURES
394 /*! Linear interpolation using exactly two FMA operations.
396 * Implements numeric equivalent of: (1-t)*d0 + t*d1
397 * Note that CUDA does not have fnms, otherwise we'd use
398 * fma(t, d1, fnms(t, d0, d0)
399 * but input modifiers are designed for this and are fast.
401 template <typename T>
402 __forceinline__ __host__ __device__
403 T lerp(T d0, T d1, T t)
405 return fma(t, d1, fma(-t, d0, d0));
408 /*! Interpolate Ewald coulomb force correction using the F*r table.
410 static __forceinline__ __device__
411 float interpolate_coulomb_force_r(const cu_nbparam_t nbparam,
414 float normalized = nbparam.coulomb_tab_scale * r;
415 int index = (int) normalized;
416 float fraction = normalized - index;
418 float2 d01 = fetch_coulomb_force_r(nbparam, index);
420 return lerp(d01.x, d01.y, fraction);
423 /*! Fetch C6 and C12 from the parameter table.
425 * Depending on what is supported, it fetches parameters either
426 * using direct load, texture objects, or texrefs.
428 static __forceinline__ __device__
429 void fetch_nbfp_c6_c12(float &c6,
431 const cu_nbparam_t nbparam,
434 #if DISABLE_CUDA_TEXTURES
435 /* Force an 8-byte fetch to save a memory instruction. */
436 float2 *nbfp = (float2 *)nbparam.nbfp;
438 c6c12 = LDG(&nbfp[baseIndex]);
442 /* NOTE: as we always do 8-byte aligned loads, we could
443 fetch float2 here too just as above. */
444 c6 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2*baseIndex);
445 c12 = tex1Dfetch<float>(nbparam.nbfp_texobj, 2*baseIndex + 1);
446 #endif // DISABLE_CUDA_TEXTURES
450 /*! Calculate analytical Ewald correction term. */
451 static __forceinline__ __device__
452 float pmecorrF(float z2)
454 const float FN6 = -1.7357322914161492954e-8f;
455 const float FN5 = 1.4703624142580877519e-6f;
456 const float FN4 = -0.000053401640219807709149f;
457 const float FN3 = 0.0010054721316683106153f;
458 const float FN2 = -0.019278317264888380590f;
459 const float FN1 = 0.069670166153766424023f;
460 const float FN0 = -0.75225204789749321333f;
462 const float FD4 = 0.0011193462567257629232f;
463 const float FD3 = 0.014866955030185295499f;
464 const float FD2 = 0.11583842382862377919f;
465 const float FD1 = 0.50736591960530292870f;
466 const float FD0 = 1.0f;
469 float polyFN0, polyFN1, polyFD0, polyFD1;
473 polyFD0 = FD4*z4 + FD2;
474 polyFD1 = FD3*z4 + FD1;
475 polyFD0 = polyFD0*z4 + FD0;
476 polyFD0 = polyFD1*z2 + polyFD0;
478 polyFD0 = 1.0f/polyFD0;
480 polyFN0 = FN6*z4 + FN4;
481 polyFN1 = FN5*z4 + FN3;
482 polyFN0 = polyFN0*z4 + FN2;
483 polyFN1 = polyFN1*z4 + FN1;
484 polyFN0 = polyFN0*z4 + FN0;
485 polyFN0 = polyFN1*z2 + polyFN0;
487 return polyFN0*polyFD0;
490 /*! Final j-force reduction; this generic implementation works with
491 * arbitrary array sizes.
493 static __forceinline__ __device__
494 void reduce_force_j_generic(float *f_buf, float3 *fout,
495 int tidxi, int tidxj, int aidx)
500 for (int j = tidxj * c_clSize; j < (tidxj + 1) * c_clSize; j++)
502 f += f_buf[c_fbufStride * tidxi + j];
505 atomicAdd((&fout[aidx].x)+tidxi, f);
509 /*! Final j-force reduction; this implementation only with power of two
510 * array sizes and with sm >= 3.0
512 #if GMX_PTX_ARCH >= 300 || GMX_PTX_ARCH == 0
513 static __forceinline__ __device__
514 void reduce_force_j_warp_shfl(float3 f, float3 *fout,
516 const unsigned int activemask)
518 f.x += gmx_shfl_down_sync(activemask, f.x, 1);
519 f.y += gmx_shfl_up_sync (activemask, f.y, 1);
520 f.z += gmx_shfl_down_sync(activemask, f.z, 1);
527 f.x += gmx_shfl_down_sync(activemask, f.x, 2);
528 f.z += gmx_shfl_up_sync (activemask, f.z, 2);
535 f.x += gmx_shfl_down_sync(activemask, f.x, 4);
539 atomicAdd((&fout[aidx].x) + tidxi, f.x);
544 /*! Final i-force reduction; this generic implementation works with
545 * arbitrary array sizes.
546 * TODO: add the tidxi < 3 trick
548 static __forceinline__ __device__
549 void reduce_force_i_generic(float *f_buf, float3 *fout,
550 float *fshift_buf, bool bCalcFshift,
551 int tidxi, int tidxj, int aidx)
556 for (int j = tidxi; j < c_clSizeSq; j += c_clSize)
558 f += f_buf[tidxj * c_fbufStride + j];
561 atomicAdd(&fout[aidx].x + tidxj, f);
570 /*! Final i-force reduction; this implementation works only with power of two
573 static __forceinline__ __device__
574 void reduce_force_i_pow2(volatile float *f_buf, float3 *fout,
575 float *fshift_buf, bool bCalcFshift,
576 int tidxi, int tidxj, int aidx)
581 assert(c_clSize == 1 << c_clSizeLog2);
583 /* Reduce the initial c_clSize values for each i atom to half
584 * every step by using c_clSize * i threads.
585 * Can't just use i as loop variable because than nvcc refuses to unroll.
589 for (j = c_clSizeLog2 - 1; j > 0; j--)
594 f_buf[ tidxj * c_clSize + tidxi] += f_buf[ (tidxj + i) * c_clSize + tidxi];
595 f_buf[ c_fbufStride + tidxj * c_clSize + tidxi] += f_buf[ c_fbufStride + (tidxj + i) * c_clSize + tidxi];
596 f_buf[2 * c_fbufStride + tidxj * c_clSize + tidxi] += f_buf[2 * c_fbufStride + (tidxj + i) * c_clSize + tidxi];
601 /* i == 1, last reduction step, writing to global mem */
604 /* tidxj*c_fbufStride selects x, y or z */
605 f = f_buf[tidxj * c_fbufStride + tidxi] +
606 f_buf[tidxj * c_fbufStride + i * c_clSize + tidxi];
608 atomicAdd(&(fout[aidx].x) + tidxj, f);
618 /*! Final i-force reduction wrapper; calls the generic or pow2 reduction depending
619 * on whether the size of the array to be reduced is power of two or not.
621 static __forceinline__ __device__
622 void reduce_force_i(float *f_buf, float3 *f,
623 float *fshift_buf, bool bCalcFshift,
624 int tidxi, int tidxj, int ai)
626 if ((c_clSize & (c_clSize - 1)))
628 reduce_force_i_generic(f_buf, f, fshift_buf, bCalcFshift, tidxi, tidxj, ai);
632 reduce_force_i_pow2(f_buf, f, fshift_buf, bCalcFshift, tidxi, tidxj, ai);
636 /*! Final i-force reduction; this implementation works only with power of two
637 * array sizes and with sm >= 3.0
639 #if GMX_PTX_ARCH >= 300 || GMX_PTX_ARCH == 0
640 static __forceinline__ __device__
641 void reduce_force_i_warp_shfl(float3 fin, float3 *fout,
642 float *fshift_buf, bool bCalcFshift,
644 const unsigned int activemask)
646 fin.x += gmx_shfl_down_sync(activemask, fin.x, c_clSize);
647 fin.y += gmx_shfl_up_sync (activemask, fin.y, c_clSize);
648 fin.z += gmx_shfl_down_sync(activemask, fin.z, c_clSize);
655 fin.x += gmx_shfl_down_sync(activemask, fin.x, 2*c_clSize);
656 fin.z += gmx_shfl_up_sync (activemask, fin.z, 2*c_clSize);
663 /* Threads 0,1,2 and 4,5,6 increment x,y,z for their warp */
666 atomicAdd(&fout[aidx].x + (tidxj & 3), fin.x);
670 *fshift_buf += fin.x;
676 /*! Energy reduction; this implementation works only with power of two
679 static __forceinline__ __device__
680 void reduce_energy_pow2(volatile float *buf,
681 float *e_lj, float *e_el,
686 unsigned int i = warp_size/2;
688 /* Can't just use i as loop variable because than nvcc refuses to unroll. */
690 for (int j = warp_size_log2 - 1; j > 0; j--)
694 buf[ tidx] += buf[ tidx + i];
695 buf[c_fbufStride + tidx] += buf[c_fbufStride + tidx + i];
700 // TODO do this on two threads - requires e_lj and e_el to be stored on adjascent
701 // memory locations to make sense
702 /* last reduction step, writing to global mem */
705 e1 = buf[ tidx] + buf[ tidx + i];
706 e2 = buf[c_fbufStride + tidx] + buf[c_fbufStride + tidx + i];
713 /*! Energy reduction; this implementation works only with power of two
714 * array sizes and with sm >= 3.0
716 #if GMX_PTX_ARCH >= 300 || GMX_PTX_ARCH == 0
717 static __forceinline__ __device__
718 void reduce_energy_warp_shfl(float E_lj, float E_el,
719 float *e_lj, float *e_el,
721 const unsigned int activemask)
727 for (i = 0; i < 5; i++)
729 E_lj += gmx_shfl_down_sync(activemask, E_lj, sh);
730 E_el += gmx_shfl_down_sync(activemask, E_el, sh);
734 /* The first thread in the warp writes the reduced energies */
735 if (tidx == 0 || tidx == warp_size)
737 atomicAdd(e_lj, E_lj);
738 atomicAdd(e_el, E_el);
741 #endif /* GMX_PTX_ARCH */
743 #endif /* NBNXN_CUDA_KERNEL_UTILS_CUH */