SYCL: Use acc.bind(cgh) instead of cgh.require(acc)
[alexxy/gromacs.git] / src / gromacs / mdlib / gpuforcereduction_impl_internal_sycl.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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  *
37  * \brief Implements GPU Force Reduction using SYCL
38  *
39  * \author Alan Gray <alang@nvidia.com>
40  * \author Andrey Alekseenko <al42and@gmail.com>
41  *
42  * \ingroup module_mdlib
43  */
44
45 #include "gmxpre.h"
46
47 #include "gpuforcereduction_impl_internal.h"
48
49 #include <utility>
50
51 #include "gromacs/gpu_utils/gmxsycl.h"
52 #include "gromacs/gpu_utils/devicebuffer.h"
53 #include "gromacs/gpu_utils/gpu_utils.h"
54 #include "gromacs/gpu_utils/gpueventsynchronizer.h"
55 #include "gromacs/utility/template_mp.h"
56
57 //! \brief Class name for reduction kernel
58 template<bool addRvecForce, bool accumulateForce>
59 class ReduceKernel;
60
61 namespace gmx
62 {
63
64 using cl::sycl::access::mode;
65
66 //! \brief Function returning the force reduction kernel lambda.
67 template<bool addRvecForce, bool accumulateForce>
68 static auto reduceKernel(cl::sycl::handler&                                 cgh,
69                          DeviceAccessor<Float3, mode::read>                 a_nbnxmForce,
70                          OptionalAccessor<Float3, mode::read, addRvecForce> a_rvecForceToAdd,
71                          DeviceAccessor<Float3, accumulateForce ? mode::read_write : mode::write> a_forceTotal,
72                          DeviceAccessor<int, cl::sycl::access::mode::read> a_cell,
73                          const int                                         atomStart)
74 {
75     a_nbnxmForce.bind(cgh);
76     if constexpr (addRvecForce)
77     {
78         a_rvecForceToAdd.bind(cgh);
79     }
80     a_forceTotal.bind(cgh);
81     a_cell.bind(cgh);
82
83     return [=](cl::sycl::id<1> itemIdx) {
84         // Set to nbnxnm force, then perhaps accumulate further to it
85         Float3 temp = a_nbnxmForce[a_cell[itemIdx]];
86
87         if constexpr (accumulateForce)
88         {
89             temp += a_forceTotal[itemIdx + atomStart];
90         }
91
92         if constexpr (addRvecForce)
93         {
94             temp += a_rvecForceToAdd[itemIdx + atomStart];
95         }
96
97         a_forceTotal[itemIdx + atomStart] = temp;
98     };
99 }
100
101 //! \brief Force reduction SYCL kernel launch code.
102 template<bool addRvecForce, bool accumulateForce>
103 static void launchReductionKernel_(const int                   numAtoms,
104                                    const int                   atomStart,
105                                    const DeviceBuffer<Float3>& b_nbnxmForce,
106                                    const DeviceBuffer<Float3>& b_rvecForceToAdd,
107                                    DeviceBuffer<Float3>&       b_forceTotal,
108                                    const DeviceBuffer<int>&    b_cell,
109                                    const DeviceStream&         deviceStream)
110 {
111     const cl::sycl::range<1> rangeNumAtoms(numAtoms);
112     cl::sycl::queue          queue = deviceStream.stream();
113
114     // We only need parts of b_rvecForceToAdd and b_forceTotal, so sub-buffers would be appropriate.
115     // But hipSYCL does not support them yet, nor plans to. See Issue #4019.
116
117     queue.submit([&](cl::sycl::handler& cgh) {
118         auto kernel = reduceKernel<addRvecForce, accumulateForce>(
119                 cgh, b_nbnxmForce, b_rvecForceToAdd, b_forceTotal, b_cell, atomStart);
120         cgh.parallel_for<ReduceKernel<addRvecForce, accumulateForce>>(rangeNumAtoms, kernel);
121     });
122 }
123
124 /*! \brief Select templated Force reduction kernel and launch it. */
125 void launchForceReductionKernel(int                  numAtoms,
126                                 int                  atomStart,
127                                 bool                 addRvecForce,
128                                 bool                 accumulate,
129                                 DeviceBuffer<Float3> d_nbnxmForceToAdd,
130                                 DeviceBuffer<Float3> d_rvecForceToAdd,
131                                 DeviceBuffer<Float3> d_baseForce,
132                                 DeviceBuffer<int>    d_cell,
133                                 const DeviceStream&  deviceStream)
134 {
135     dispatchTemplatedFunction(
136             [&](auto addRvecForce_, auto accumulateForce_) {
137                 return launchReductionKernel_<addRvecForce_, accumulateForce_>(
138                         numAtoms, atomStart, d_nbnxmForceToAdd, d_rvecForceToAdd, d_baseForce, d_cell, deviceStream);
139             },
140             addRvecForce,
141             accumulate);
142 }
143
144 } // namespace gmx