c3aae4ebb50d930b7e14cbb506ba04ac91ba5a8c
[alexxy/gromacs.git] / src / gromacs / nbnxm / sycl / nbnxm_sycl_kernel_utils.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020,2021, 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 /*! \internal \file
36  * \brief
37  * Helper functions and constants for SYCL NBNXM kernels
38  *
39  * \ingroup module_nbnxm
40  */
41 #ifndef GMX_NBNXM_SYCL_NBNXN_SYCL_KERNEL_UTILS_H
42 #define GMX_NBNXM_SYCL_NBNXN_SYCL_KERNEL_UTILS_H
43
44 #include "gromacs/nbnxm/pairlist.h"
45 #include "gromacs/nbnxm/pairlistparams.h"
46
47 namespace Nbnxm
48 {
49
50 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
51 #    define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY 4
52 #endif
53 /*! \brief Macro defining default for the prune kernel's j4 processing concurrency.
54  *
55  *  The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
56  */
57 static constexpr int c_syclPruneKernelJ4Concurrency = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY;
58
59 /* Convenience constants */
60 /*! \cond */
61 // cluster size = number of atoms per cluster.
62 static constexpr int c_clSize = c_nbnxnGpuClusterSize;
63 // j-cluster size after split (4 in the current implementation).
64 static constexpr int c_splitClSize = c_clSize / c_nbnxnGpuClusterpairSplit;
65 // i-cluster interaction mask for a super-cluster with all c_nbnxnGpuNumClusterPerSupercluster=8 bits set.
66 static constexpr unsigned superClInteractionMask = ((1U << c_nbnxnGpuNumClusterPerSupercluster) - 1U);
67
68 // 1/sqrt(pi), same value as \c M_FLOAT_1_SQRTPI in other NB kernels.
69 static constexpr float c_OneOverSqrtPi = 0.564189583547756F;
70 // 1/6, same value as in other NB kernels.
71 static constexpr float c_oneSixth = 0.16666667F;
72 // 1/12, same value as in other NB kernels.
73 static constexpr float c_oneTwelfth = 0.08333333F;
74 /*! \endcond */
75
76 /* The following functions are necessary because on some versions of Intel OpenCL RT, subgroups
77  * do not properly work (segfault or create subgroups of size 1) if used in kernels
78  * with non-1-dimensional workgroup. */
79 //! \brief Convert 3D range to 1D
80 static inline cl::sycl::range<1> flattenRange(cl::sycl::range<3> range3d)
81 {
82     return cl::sycl::range<1>(range3d.size());
83 }
84
85 //! \brief Convert 3D nd_range to 1D
86 static inline cl::sycl::nd_range<1> flattenNDRange(cl::sycl::nd_range<3> nd_range3d)
87 {
88     return cl::sycl::nd_range<1>(flattenRange(nd_range3d.get_global_range()),
89                                  flattenRange(nd_range3d.get_local_range()));
90 }
91
92 //! \brief Convert flattened 1D index to 3D
93 template<int rangeX, int rangeY>
94 static inline cl::sycl::id<3> unflattenId(cl::sycl::id<1> id1d)
95 {
96     constexpr unsigned rangeXY = rangeX * rangeY;
97     const unsigned     id      = id1d[0];
98     const unsigned     z       = id / rangeXY;
99     const unsigned     xy      = id % rangeXY;
100     return cl::sycl::id<3>(xy % rangeX, xy / rangeX, z);
101 }
102
103 //! \brief Convenience wrapper to do atomic addition to a global buffer
104 template<cl::sycl::access::mode Mode, class IndexType>
105 static inline void atomicFetchAdd(DeviceAccessor<float, Mode> acc, const IndexType idx, const float val)
106 {
107     if (cl::sycl::isnormal(val))
108     {
109         sycl_2020::atomic_ref<float, sycl_2020::memory_order::relaxed, sycl_2020::memory_scope::device, cl::sycl::access::address_space::global_space>
110                 fout_atomic(acc[idx]);
111         fout_atomic.fetch_add(val);
112     }
113 }
114
115 static inline float shuffleDown(float var, unsigned int delta, sycl_2020::sub_group sg)
116 {
117     return sg.shuffle_down(var, delta);
118 }
119
120 static inline float shuffleUp(float var, unsigned int delta, sycl_2020::sub_group sg)
121 {
122     return sg.shuffle_up(var, delta);
123 }
124
125 } // namespace Nbnxm
126
127 #endif // GMX_NBNXM_SYCL_NBNXN_SYCL_KERNEL_UTILS_H