Add common header to GpuEventSynchronizer
[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 template<bool addRvecForce, bool accumulateForce>
67 static auto reduceKernel(cl::sycl::handler&                                 cgh,
68                          DeviceAccessor<Float3, mode::read>                 a_nbnxmForce,
69                          OptionalAccessor<Float3, mode::read, addRvecForce> a_rvecForceToAdd,
70                          DeviceAccessor<Float3, accumulateForce ? mode::read_write : mode::write> a_forceTotal,
71                          DeviceAccessor<int, cl::sycl::access::mode::read> a_cell,
72                          const int                                         atomStart)
73 {
74     cgh.require(a_nbnxmForce);
75     if constexpr (addRvecForce)
76     {
77         cgh.require(a_rvecForceToAdd);
78     }
79     cgh.require(a_forceTotal);
80     cgh.require(a_cell);
81
82     return [=](cl::sycl::id<1> itemIdx) {
83         // Set to nbnxnm force, then perhaps accumulate further to it
84         Float3 temp = a_nbnxmForce[a_cell[itemIdx]];
85
86         if constexpr (accumulateForce)
87         {
88             temp += a_forceTotal[itemIdx + atomStart];
89         }
90
91         if constexpr (addRvecForce)
92         {
93             temp += a_rvecForceToAdd[itemIdx + atomStart];
94         }
95
96         a_forceTotal[itemIdx + atomStart] = temp;
97     };
98 }
99
100 template<bool addRvecForce, bool accumulateForce>
101 static void launchReductionKernel_(const int                   numAtoms,
102                                    const int                   atomStart,
103                                    const DeviceBuffer<Float3>& b_nbnxmForce,
104                                    const DeviceBuffer<Float3>& b_rvecForceToAdd,
105                                    DeviceBuffer<Float3>&       b_forceTotal,
106                                    const DeviceBuffer<int>&    b_cell,
107                                    const DeviceStream&         deviceStream)
108 {
109     const cl::sycl::range<1> rangeNumAtoms(numAtoms);
110     cl::sycl::queue          queue = deviceStream.stream();
111
112     // We only need parts of b_rvecForceToAdd and b_forceTotal, so sub-buffers would be appropriate.
113     // But hipSYCL does not support them yet, nor plans to. See Issue #4019.
114
115     queue.submit([&](cl::sycl::handler& cgh) {
116         auto kernel = reduceKernel<addRvecForce, accumulateForce>(
117                 cgh, b_nbnxmForce, b_rvecForceToAdd, b_forceTotal, b_cell, atomStart);
118         cgh.parallel_for<ReduceKernel<addRvecForce, accumulateForce>>(rangeNumAtoms, kernel);
119     });
120 }
121
122 /*! \brief Select templated kernel and launch it. */
123 void launchForceReductionKernel(int                  numAtoms,
124                                 int                  atomStart,
125                                 bool                 addRvecForce,
126                                 bool                 accumulate,
127                                 DeviceBuffer<Float3> d_nbnxmForceToAdd,
128                                 DeviceBuffer<Float3> d_rvecForceToAdd,
129                                 DeviceBuffer<Float3> d_baseForce,
130                                 DeviceBuffer<int>    d_cell,
131                                 const DeviceStream&  deviceStream)
132 {
133     dispatchTemplatedFunction(
134             [&](auto addRvecForce_, auto accumulateForce_) {
135                 return launchReductionKernel_<addRvecForce_, accumulateForce_>(
136                         numAtoms, atomStart, d_nbnxmForceToAdd, d_rvecForceToAdd, d_baseForce, d_cell, deviceStream);
137             },
138             addRvecForce,
139             accumulate);
140 }
141
142 } // namespace gmx