Add SYCL implementation of GPU F buffer operations
authorAndrey Alekseenko <al42and@gmail.com>
Tue, 11 May 2021 07:16:38 +0000 (07:16 +0000)
committerArtem Zhmurov <zhmurov@gmail.com>
Tue, 11 May 2021 07:16:38 +0000 (07:16 +0000)
src/gromacs/mdlib/CMakeLists.txt
src/gromacs/mdlib/gpuforcereduction.h
src/gromacs/mdlib/gpuforcereduction_impl.cpp
src/gromacs/mdlib/gpuforcereduction_impl_internal_sycl.cpp [new file with mode: 0644]

index 060f6160d7ba9dea6d9de1bdb53eb6beb001b449..511a94680966be1e93ce89ed6a5b9cefe62d9aef 100644 (file)
@@ -39,6 +39,7 @@ file(GLOB MDLIB_SOURCES *.cpp)
 # To avoid listing all the necessary files manually, we will remove SYCL-specific files here:
 list(REMOVE_ITEM MDLIB_SOURCES
     ${CMAKE_CURRENT_SOURCE_DIR}/gpuforcereduction_impl.cpp
+    ${CMAKE_CURRENT_SOURCE_DIR}/gpuforcereduction_impl_internal_sycl.cpp
     ${CMAKE_CURRENT_SOURCE_DIR}/leapfrog_gpu_sycl.cpp
     ${CMAKE_CURRENT_SOURCE_DIR}/lincs_gpu.cpp
     ${CMAKE_CURRENT_SOURCE_DIR}/lincs_gpu_internal_sycl.cpp
@@ -66,6 +67,8 @@ endif()
 
 if(GMX_GPU_SYCL)
     gmx_add_libgromacs_sources(
+        gpuforcereduction_impl.cpp
+        gpuforcereduction_impl_internal_sycl.cpp
         leapfrog_gpu_sycl.cpp
         lincs_gpu.cpp
         lincs_gpu_internal_sycl.cpp
@@ -74,6 +77,8 @@ if(GMX_GPU_SYCL)
     )
 
     _gmx_add_files_to_property(SYCL_SOURCES
+        gpuforcereduction_impl.cpp
+        gpuforcereduction_impl_internal_sycl.cpp
         leapfrog_gpu_sycl.cpp
         lincs_gpu.cpp
         lincs_gpu_internal_sycl.cpp
index e015d7ef5e03c70b4782ee74f7e8b13e622bbaea..922b8fecaf4fd7a2dd432fd1e9cf02f37c810f44 100644 (file)
@@ -60,7 +60,7 @@ class DeviceContext;
 namespace gmx
 {
 
-#define HAVE_GPU_FORCE_REDUCTION (GMX_GPU_CUDA)
+#define HAVE_GPU_FORCE_REDUCTION (GMX_GPU_CUDA || GMX_GPU_SYCL)
 
 /*! \internal
  * \brief Manages the force reduction directly in GPU memory
index fb58c5c9437ef48cbadfb5ac4c9592d9b4f006ea..93772853d5301781371172fe4157cb86514a844e 100644 (file)
@@ -77,7 +77,7 @@ void GpuForceReduction::Impl::reinit(DeviceBuffer<Float3>  baseForcePtr,
                                      const bool            accumulate,
                                      GpuEventSynchronizer* completionMarker)
 {
-    GMX_ASSERT((baseForcePtr != nullptr), "Input base force for reduction has no data");
+    GMX_ASSERT(baseForcePtr, "Input base force for reduction has no data");
     baseForce_        = baseForcePtr;
     numAtoms_         = numAtoms;
     atomStart_        = atomStart;
diff --git a/src/gromacs/mdlib/gpuforcereduction_impl_internal_sycl.cpp b/src/gromacs/mdlib/gpuforcereduction_impl_internal_sycl.cpp
new file mode 100644 (file)
index 0000000..3bc21c4
--- /dev/null
@@ -0,0 +1,141 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2021, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ *
+ * \brief Implements GPU Force Reduction using SYCL
+ *
+ * \author Alan Gray <alang@nvidia.com>
+ * \author Andrey Alekseenko <al42and@gmail.com>
+ *
+ * \ingroup module_mdlib
+ */
+
+#include "gmxpre.h"
+
+#include "gpuforcereduction_impl_internal.h"
+
+#include <utility>
+
+#include "gromacs/gpu_utils/gmxsycl.h"
+#include "gromacs/gpu_utils/devicebuffer.h"
+#include "gromacs/gpu_utils/gpu_utils.h"
+#include "gromacs/gpu_utils/gpueventsynchronizer_sycl.h"
+#include "gromacs/utility/template_mp.h"
+
+namespace gmx
+{
+
+using cl::sycl::access::mode;
+
+template<bool addRvecForce, bool accumulateForce>
+static auto reduceKernel(cl::sycl::handler&                                 cgh,
+                         DeviceAccessor<Float3, mode::read>                 a_nbnxmForce,
+                         OptionalAccessor<Float3, mode::read, addRvecForce> a_rvecForceToAdd,
+                         DeviceAccessor<Float3, accumulateForce ? mode::read_write : mode::write> a_forceTotal,
+                         DeviceAccessor<int, cl::sycl::access::mode::read> a_cell,
+                         const int                                         atomStart)
+{
+    cgh.require(a_nbnxmForce);
+    if constexpr (addRvecForce)
+    {
+        cgh.require(a_rvecForceToAdd);
+    }
+    cgh.require(a_forceTotal);
+    cgh.require(a_cell);
+
+    return [=](cl::sycl::id<1> itemIdx) {
+        // Set to nbnxnm force, then perhaps accumulate further to it
+        Float3 temp = a_nbnxmForce[a_cell[itemIdx]];
+
+        if constexpr (accumulateForce)
+        {
+            temp += a_forceTotal[itemIdx + atomStart];
+        }
+
+        if constexpr (addRvecForce)
+        {
+            temp += a_rvecForceToAdd[itemIdx + atomStart];
+        }
+
+        a_forceTotal[itemIdx + atomStart] = temp;
+    };
+}
+
+template<bool addRvecForce, bool accumulateForce>
+class ReduceKernelName;
+
+template<bool addRvecForce, bool accumulateForce>
+static void launchReductionKernel_(const int                   numAtoms,
+                                   const int                   atomStart,
+                                   const DeviceBuffer<Float3>& b_nbnxmForce,
+                                   const DeviceBuffer<Float3>& b_rvecForceToAdd,
+                                   DeviceBuffer<Float3>&       b_forceTotal,
+                                   const DeviceBuffer<int>&    b_cell,
+                                   const DeviceStream&         deviceStream)
+{
+    const cl::sycl::range<1> rangeNumAtoms(numAtoms);
+    cl::sycl::queue          queue = deviceStream.stream();
+
+    // We only need parts of b_rvecForceToAdd and b_forceTotal, so sub-buffers would be appropriate.
+    // But hipSYCL does not support them yet, nor plans to. See Issue #4019.
+
+    queue.submit([&](cl::sycl::handler& cgh) {
+        auto kernel = reduceKernel<addRvecForce, accumulateForce>(
+                cgh, b_nbnxmForce, b_rvecForceToAdd, b_forceTotal, b_cell, atomStart);
+        cgh.parallel_for<ReduceKernelName<addRvecForce, accumulateForce>>(rangeNumAtoms, kernel);
+    });
+}
+
+/*! \brief Select templated kernel and launch it. */
+void launchForceReductionKernel(int                  numAtoms,
+                                int                  atomStart,
+                                bool                 addRvecForce,
+                                bool                 accumulate,
+                                DeviceBuffer<Float3> d_nbnxmForceToAdd,
+                                DeviceBuffer<Float3> d_rvecForceToAdd,
+                                DeviceBuffer<Float3> d_baseForce,
+                                DeviceBuffer<int>    d_cell,
+                                const DeviceStream&  deviceStream)
+{
+    dispatchTemplatedFunction(
+            [&](auto addRvecForce_, auto accumulateForce_) {
+                return launchReductionKernel_<addRvecForce_, accumulateForce_>(
+                        numAtoms, atomStart, d_nbnxmForceToAdd, d_rvecForceToAdd, d_baseForce, d_cell, deviceStream);
+            },
+            addRvecForce,
+            accumulate);
+}
+
+} // namespace gmx