2 * This file is part of the GROMACS molecular simulation package.
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.
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.
37 * \brief Implements GPU Force Reduction using SYCL
39 * \author Alan Gray <alang@nvidia.com>
40 * \author Andrey Alekseenko <al42and@gmail.com>
42 * \ingroup module_mdlib
47 #include "gpuforcereduction_impl_internal.h"
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"
60 using cl::sycl::access::mode;
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,
70 cgh.require(a_nbnxmForce);
71 if constexpr (addRvecForce)
73 cgh.require(a_rvecForceToAdd);
75 cgh.require(a_forceTotal);
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]];
82 if constexpr (accumulateForce)
84 temp += a_forceTotal[itemIdx + atomStart];
87 if constexpr (addRvecForce)
89 temp += a_rvecForceToAdd[itemIdx + atomStart];
92 a_forceTotal[itemIdx + atomStart] = temp;
96 template<bool addRvecForce, bool accumulateForce>
97 class ReduceKernelName;
99 template<bool addRvecForce, bool accumulateForce>
100 static void launchReductionKernel_(const int numAtoms,
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)
108 const cl::sycl::range<1> rangeNumAtoms(numAtoms);
109 cl::sycl::queue queue = deviceStream.stream();
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.
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);
121 /*! \brief Select templated kernel and launch it. */
122 void launchForceReductionKernel(int numAtoms,
126 DeviceBuffer<Float3> d_nbnxmForceToAdd,
127 DeviceBuffer<Float3> d_rvecForceToAdd,
128 DeviceBuffer<Float3> d_baseForce,
129 DeviceBuffer<int> d_cell,
130 const DeviceStream& deviceStream)
132 dispatchTemplatedFunction(
133 [&](auto addRvecForce_, auto accumulateForce_) {
134 return launchReductionKernel_<addRvecForce_, accumulateForce_>(
135 numAtoms, atomStart, d_nbnxmForceToAdd, d_rvecForceToAdd, d_baseForce, d_cell, deviceStream);