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