Enable GMX_GPU=SYCL compiles
authorErik Lindahl <erik.lindahl@gmail.com>
Tue, 15 Sep 2020 12:38:19 +0000 (12:38 +0000)
committerPaul Bauer <paul.bauer.q@gmail.com>
Tue, 15 Sep 2020 12:38:19 +0000 (12:38 +0000)
Update CMake config to properly apply -fsycl, but only to files
including actual SYCL code, and get linking working.  Tested
with both dpcpp and icpx, and added trivial SYCL code to a
dimmy device context, just to ensure compiles and linking works
both with and without SYCL support enabled.

Also fixes a number of new warnings about unused variables now
detected with dpcpp and icpx.

CMakeLists.txt
cmake/gmxCFlags.cmake
cmake/gmxManageSYCL.cmake
src/gromacs/CMakeLists.txt
src/gromacs/gpu_utils/CMakeLists.txt
src/gromacs/gpu_utils/clfftinitializer.cpp
src/gromacs/gpu_utils/device_context_sycl.cpp [new file with mode: 0644]
src/gromacs/gpu_utils/device_stream_sycl.cpp [new file with mode: 0644]
src/gromacs/gpu_utils/devicebuffer_sycl.h
src/gromacs/gpu_utils/gpueventsynchronizer_sycl.h

index 58f79778abceebf05ba6d988e40a52fde979546e..646042333e5a43c7a821b7e23ddcce971a6c64e3 100644 (file)
@@ -716,15 +716,15 @@ endif()
 # # # # # # # # # # NO MORE TESTS AFTER THIS LINE! # # # # # # # # # # #
 # these are set after everything else
 if (NOT GMX_SKIP_DEFAULT_CFLAGS)
-    set(CMAKE_EXE_LINKER_FLAGS "${FFT_LINKER_FLAGS} ${MPI_LINKER_FLAGS} ${CMAKE_EXE_LINKER_FLAGS}")
-    set(CMAKE_SHARED_LINKER_FLAGS "${FFT_LINKER_FLAGS} ${MPI_LINKER_FLAGS} ${CMAKE_SHARED_LINKER_FLAGS}")
+    set(CMAKE_EXE_LINKER_FLAGS "${FFT_LINKER_FLAGS} ${MPI_LINKER_FLAGS} ${CMAKE_EXE_LINKER_FLAGS} ${DISABLE_SYCL_CXX_FLAGS}")
+    set(CMAKE_SHARED_LINKER_FLAGS "${FFT_LINKER_FLAGS} ${MPI_LINKER_FLAGS} ${CMAKE_SHARED_LINKER_FLAGS} ${DISABLE_SYCL_CXX_FLAGS}")
 else()
     message("Recommended flags which are not added because GMX_SKIP_DEFAULT_CFLAGS=yes:")
     message("CMAKE_C_FLAGS: ${SIMD_C_FLAGS};${MPI_COMPILE_FLAGS};${EXTRA_C_FLAGS};${GMXC_CFLAGS}")
     foreach(build_type ${build_types_with_explicit_flags})
         message("CMAKE_C_FLAGS_${build_type}: ${GMXC_CFLAGS_${build_type}}")
     endforeach()
-    message("CMAKE_CXX_FLAGS: ${DISABLE_SYCL_CXX_FLAGS};${SIMD_CXX_FLAGS};${MPI_COMPILE_FLAGS};${EXTRA_CXX_FLAGS};${GMXC_CXXFLAGS}")
+    message("CMAKE_CXX_FLAGS: ${SIMD_CXX_FLAGS};${MPI_COMPILE_FLAGS};${EXTRA_CXX_FLAGS};${GMXC_CXXFLAGS}")
     foreach(build_type ${build_types_with_explicit_flags})
         message("CMAKE_CXX_FLAGS_${build_type}: ${GMXC_CXXFLAGS_${build_type}}")
     endforeach()
index 5789f62e51cdfb884d4d3932e85b3873f883078b..1fffd3e7a90f2228adbf0a3d82faad8adeba7dfe 100644 (file)
@@ -60,7 +60,10 @@ ENDMACRO(GMX_TEST_CXXFLAG VARIABLE FLAGS CXXFLAGSVAR)
 # works the same way.
 function(gmx_target_compile_options_inner)
     set (CFLAGS "${SIMD_C_FLAGS};${MPI_COMPILE_FLAGS};${EXTRA_C_FLAGS};${GMXC_CFLAGS}" PARENT_SCOPE)
-    set (CXXFLAGS "${SIMD_CXX_FLAGS};${MPI_COMPILE_FLAGS};${EXTRA_CXX_FLAGS};${GMXC_CXXFLAGS}" PARENT_SCOPE)
+    # When SYCL support has been enabled (so the flag is non-empty), we still *disable* things
+    # by default to avoid running each file three passes through the compiler. Then we'll explicitly
+    # enable SYCL for the few files using it, as well as the linker.
+    set (CXXFLAGS "${SIMD_CXX_FLAGS};${MPI_COMPILE_FLAGS};${DISABLE_SYCL_CXX_FLAGS};${EXTRA_CXX_FLAGS};${GMXC_CXXFLAGS}" PARENT_SCOPE)
 endfunction()
 
 # Implementation function to add compiler flags expected for all
index 0cbc8a5ec96935d1f3b3ad2d917a1084de1531d8..e19cafb6842f41d2111aa515fa293ffb6b987b77 100644 (file)
@@ -57,22 +57,35 @@ if(CMAKE_CXX_COMPILER MATCHES "dpcpp")
     # a failed link stage (when the same flag is not used). Since none of this is critical, we handle
     # it by merely checking if it works to compile a source fils with this flag, and choking if SYCL
     # is still enabled.
-    message(STATUS "Checking for flags to disable SYCL for non-SYCL hardware-specific C++ code")
+
+    if(NOT CHECK_DISABLE_SYCL_CXX_FLAGS_QUIETLY)
+        message(STATUS "Checking for flags to disable SYCL")
+    endif()
+
     gmx_check_source_compiles_with_flags(
         "int main() { return 0; }"
         "-fno-sycl"
        "CXX"
-        DISABLE_SYCL_CXX_FLAG_RESULT)
-    if(DISABLE_SYCL_CXX_FLAG_RESULT)
-        set(DISABLE_SYCL_FLAGS "-fno-sycl")
-    else()
-        message(WARNING "Cannot find flags to disable SYCL for non-SYCL hardware-specific C++ code. Expect many warnings, but they are likely benign.")
+        DISABLE_SYCL_CXX_FLAGS_RESULT)
+
+    if(DISABLE_SYCL_CXX_FLAGS_RESULT)
+        set(DISABLE_SYCL_CXX_FLAGS "-fno-sycl")
+    endif()
+    if(NOT CHECK_DISABLE_SYCL_CXX_FLAGS_QUIETLY)
+        if(DISABLE_SYCL_CXX_FLAGS_RESULT)
+            message(STATUS "Checking for flags to disable SYCL - -fno-sycl")
+        else()
+            message(WARNING "Cannot find flags to disable SYCL for non-SYCL hardware-specific C++ code. Expect many warnings, but they are likely benign.")
+        endif()
+        set(CHECK_DISABLE_SYCL_CXX_FLAGS_QUIETLY 1 CACHE INTERNAL "Keep quiet on future calls to detect no-SYCL flags" FORCE)
     endif()
 endif()
+
 # Find the flags to enable (or re-enable) SYCL with Intel extensions. In case we turned it off above,
 # it's important that we check the combination of both flags, to make sure the second one re-enables SYCL.
-message(STATUS "Checking for flags to enable SYCL source")
+if(NOT CHECK_SYCL_CXX_FLAGS_QUIETLY)
+    message(STATUS "Checking for flags to enable SYCL")
+endif()
 gmx_find_flag_for_source(SYCL_CXX_FLAGS_RESULT
     "#include <CL/sycl.hpp>
      namespace sycl = cl::sycl;
@@ -87,10 +100,15 @@ gmx_find_flag_for_source(SYCL_CXX_FLAGS_RESULT
          q.parallel_for<class whatever>(sycl::range<1>{length}, [=, ptr = v.data()] (sycl::id<1> i){ ptr[i] *= 2; }).wait();
          return 0;
      }
-     " "CXX" DISABLE_SYCL_FLAGS SYCL_CXX_FLAGS "-fsycl")
-  
-if(SYCL_CXX_FLAGS_RESULT)
-    message(STATUS "Checking for flags to enable SYCL source - found")
-else()
+     " "CXX" DISABLE_SYCL_CXX_FLAGS SYCL_CXX_FLAGS "-fsycl")
+
+if(NOT CHECK_SYCL_CXX_FLAGS_QUIETLY)
+    if(SYCL_CXX_FLAGS_RESULT)
+        message(STATUS "Checking for flags to enable SYCL - ${SYCL_CXX_FLAGS}")
+    endif()
+    set(CHECK_SYCL_CXX_FLAGS_QUIETLY 1 CACHE INTERNAL "Keep quiet on future calls to detect SYCL flags" FORCE)
+endif()
+
+if(NOT SYCL_CXX_FLAGS_RESULT)
     message(ERROR "Cannot compile SYCL with Intel extensions. Try a different compiler or disable SYCL.")
 endif()
index 49bf5c59d9ac896df0373e5eb0ad8d6148b99cb7..ec0422488435cd7115fa6f65af02e70abd36b985 100644 (file)
@@ -282,8 +282,16 @@ if(SIMD_AVX_512_CXX_SUPPORTED AND NOT ("${GMX_SIMD_ACTIVE}" STREQUAL "AVX_512_KN
     set_source_files_properties(hardware/identifyavx512fmaunits.cpp PROPERTIES COMPILE_FLAGS "${SIMD_AVX_512_CXX_FLAGS} ${CXX_NO_UNUSED_OPTION_WARNING_FLAGS}")
 endif()
 
+# Only add the -fsycl flag to sources that really need it
+get_property(SYCL_SOURCES GLOBAL PROPERTY SYCL_SOURCES)
+set_source_files_properties(${SYCL_SOURCES} PROPERTIES COMPILE_FLAGS "${SYCL_CXX_FLAGS}")
+
 gmx_setup_tng_for_libgromacs()
 
+# We apply the SYCL flag explicitly just for libgromacs, since bugs in the beta versions of
+# icpx/dpcpp leads to crashes if we try to link an library without any SYCL code with the
+# -fsycl flag enabled. Once that bug is fixed, we should change it to simply add
+# SYCL_CXX_FLAGS to GMX_SHARED_LINKER_FLAGS.
 target_link_libraries(libgromacs
                       PRIVATE
                       ${EXTRAE_LIBRARIES}
@@ -291,6 +299,7 @@ target_link_libraries(libgromacs
                       ${GMX_COMMON_LIBRARIES}
                       ${FFT_LIBRARIES} ${LINEAR_ALGEBRA_LIBRARIES}
                       ${THREAD_LIB} ${GMX_SHARED_LINKER_FLAGS}
+                     ${SYCL_CXX_FLAGS}
                       ${OpenCL_LIBRARIES}
                       $<$<PLATFORM_ID:SunOS>:socket>
                       PUBLIC
index 58a81157d821c73599b771a082f7f4c0c0163284..a9b1a89f52f494687a2f29026550ab13a589b344 100644 (file)
@@ -58,6 +58,12 @@ elseif(GMX_GPU_CUDA)
         pinning.cu
         pmalloc_cuda.cu
         )
+elseif(GMX_GPU_SYCL)
+    gmx_add_libgromacs_sources(
+        device_context_sycl.cpp
+        device_stream_sycl.cpp
+        )
+    _gmx_add_files_to_property(SYCL_SOURCES device_context_sycl.cpp device_stream_sycl.cpp)
 else()
     gmx_add_libgromacs_sources(
         device_context.cpp
index bcaf26d56dde35f5309577718a38f95ee8f80e60..5097649d478cfd8608375da038479712e52f8cb2 100644 (file)
@@ -67,8 +67,8 @@ namespace
  * This ensures that thread-MPI and OpenMP builds can't accidentally
  * initialize it more than once. */
 //! @{
-gmx_unused bool g_clfftInitialized = false;
-gmx::Mutex      g_clfftMutex;
+bool       g_clfftInitialized = false;
+gmx::Mutex g_clfftMutex;
 //! @}
 #endif
 
diff --git a/src/gromacs/gpu_utils/device_context_sycl.cpp b/src/gromacs/gpu_utils/device_context_sycl.cpp
new file mode 100644 (file)
index 0000000..06f78f7
--- /dev/null
@@ -0,0 +1,70 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2020, 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 the DeviceContext for SYCL builds.
+ *
+ * \author Erik Lindahl <erik.lindahl@gmail.com>
+ *
+ * \ingroup module_gpu_utils
+ */
+#include "gmxpre.h"
+
+#include <vector>
+
+#include <CL/sycl.hpp>
+
+#include "gromacs/gpu_utils/device_context.h"
+
+
+//! Constructor.
+DeviceContext::DeviceContext(const DeviceInformation& deviceInfo) : deviceInfo_(deviceInfo)
+{
+    // This code is just a meaningless placeholder for now, but to actually test that
+    // compilation with and without SYCL works, and that the correct flags are set, it's
+    // important to have a file that actually contains SYCL code.
+    std::vector<cl::sycl::platform> platforms = cl::sycl::platform::get_platforms();
+    if (platforms.size() > 0)
+    {
+        for (const auto& platform : platforms)
+        {
+            printf("Found SYCL platform %s.\n",
+                   platform.get_info<cl::sycl::info::platform::name>().c_str());
+        }
+    }
+}
+
+//! Destructor
+DeviceContext::~DeviceContext() = default;
diff --git a/src/gromacs/gpu_utils/device_stream_sycl.cpp b/src/gromacs/gpu_utils/device_stream_sycl.cpp
new file mode 100644 (file)
index 0000000..275342e
--- /dev/null
@@ -0,0 +1,63 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2020, 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 the DeviceStream for SYCL builds.
+ *
+ * \author Erik Lindahl <erik.lindahl@gmail.com>
+ *
+ * \ingroup module_gpu_utils
+ */
+#include "gmxpre.h"
+
+#include "gromacs/gpu_utils/device_stream.h"
+
+DeviceStream::DeviceStream() = default;
+
+void DeviceStream::init(const DeviceContext& /* deviceContext */,
+                        DeviceStreamPriority /* priority */,
+                        const bool /* useTiming */)
+{
+}
+
+DeviceStream::~DeviceStream() = default;
+
+// NOLINTNEXTLINE readability-convert-member-functions-to-static
+bool DeviceStream::isValid() const
+{
+    return false;
+}
+
+void DeviceStream::synchronize() const {};
index bb2bc118fc2adaac22b82b28a0bdac7ba633a6ad..e487160b232e7ab0a87d48accdfe6f9aca96a0ee 100644 (file)
  * It is currently a caller's responsibility to call it only on not-yet allocated buffers.
  *
  * \tparam        ValueType            Raw value type of the \p buffer.
- * param[in,out] buffer               Pointer to the device-side buffer.
- * param[in]     numValues            Number of values to accomodate.
- * param[in]     deviceContext        The buffer's device context-to-be.
+ * \param[in,out] buffer               Pointer to the device-side buffer.
+ * \param[in]     numValues            Number of values to accomodate.
+ * \param[in]     deviceContext        The buffer's device context-to-be.
  */
 template<typename ValueType>
-void allocateDeviceBuffer(DeviceBuffer<ValueType>* /* buffer */,
-                          size_t /* numValues */,
-                          const DeviceContext& /* deviceContext */)
+void allocateDeviceBuffer(DeviceBuffer<ValueType>* gmx_unused buffer,
+                          size_t gmx_unused    numValues,
+                          const DeviceContext& gmx_unused deviceContext)
 {
     // SYCL-TODO
 }
@@ -76,10 +76,10 @@ void allocateDeviceBuffer(DeviceBuffer<ValueType>* /* buffer */,
  * as this is planned to be a destructor of DeviceBuffer as a proper class,
  * and no calls on \p buffer should be made afterwards.
  *
- * param[in] buffer  Pointer to the buffer to free.
+ * \param[in] buffer  Pointer to the buffer to free.
  */
 template<typename DeviceBuffer>
-void freeDeviceBuffer(DeviceBuffer* /* buffer */)
+void freeDeviceBuffer(DeviceBuffer* gmx_unused buffer)
 {
     // SYCL-TODO
 }
@@ -91,24 +91,24 @@ void freeDeviceBuffer(DeviceBuffer* /* buffer */)
  * because of the early return.
  *
  * \tparam        ValueType            Raw value type of the \p buffer.
- * param[in,out] buffer               Pointer to the device-side buffer
- * param[in]     hostBuffer           Pointer to the raw host-side memory, also typed \p ValueType
- * param[in]     startingOffset       Offset (in values) at the device-side buffer to copy into.
- * param[in]     numValues            Number of values to copy.
- * param[in]     deviceStream         GPU stream to perform asynchronous copy in.
- * param[in]     transferKind         Copy type: synchronous or asynchronous.
- * param[out]    timingEvent          A pointer to the H2D copy timing event to be filled in.
+ * \param[in,out] buffer               Pointer to the device-side buffer
+ * \param[in]     hostBuffer           Pointer to the raw host-side memory, also typed \p ValueType
+ * \param[in]     startingOffset       Offset (in values) at the device-side buffer to copy into.
+ * \param[in]     numValues            Number of values to copy.
+ * \param[in]     deviceStream         GPU stream to perform asynchronous copy in.
+ * \param[in]     transferKind         Copy type: synchronous or asynchronous.
+ * \param[out]    timingEvent          A pointer to the H2D copy timing event to be filled in.
  *                                     If the pointer is not null, the event can further be used
  *                                     to queue a wait for this operation or to query profiling information.
  */
 template<typename ValueType>
-void copyToDeviceBuffer(DeviceBuffer<ValueType>* /* buffer */,
-                        const ValueType* /* hostBuffer */,
-                        size_t /* startingOffset */,
-                        size_t /* numValues */,
-                        const DeviceStream& /* deviceStream */,
-                        GpuApiCallBehavior /* transferKind */,
-                        CommandEvent* /* timingEvent */)
+void copyToDeviceBuffer(DeviceBuffer<ValueType>* gmx_unused buffer,
+                        const ValueType* gmx_unused hostBuffer,
+                        size_t gmx_unused startingOffset,
+                        size_t gmx_unused   numValues,
+                        const DeviceStream& gmx_unused deviceStream,
+                        GpuApiCallBehavior gmx_unused transferKind,
+                        CommandEvent* gmx_unused timingEvent)
 {
     // SYCL-TODO
 }
@@ -120,24 +120,24 @@ void copyToDeviceBuffer(DeviceBuffer<ValueType>* /* buffer */,
  * because of the early return.
  *
  * \tparam        ValueType            Raw value type of the \p buffer.
- * param[in,out] hostBuffer           Pointer to the raw host-side memory, also typed \p ValueType
- * param[in]     buffer               Pointer to the device-side buffer
- * param[in]     startingOffset       Offset (in values) at the device-side buffer to copy from.
- * param[in]     numValues            Number of values to copy.
- * param[in]     deviceStream         GPU stream to perform asynchronous copy in.
- * param[in]     transferKind         Copy type: synchronous or asynchronous.
- * param[out]    timingEvent          A pointer to the H2D copy timing event to be filled in.
+ * \param[in,out] hostBuffer           Pointer to the raw host-side memory, also typed \p ValueType
+ * \param[in]     buffer               Pointer to the device-side buffer
+ * \param[in]     startingOffset       Offset (in values) at the device-side buffer to copy from.
+ * \param[in]     numValues            Number of values to copy.
+ * \param[in]     deviceStream         GPU stream to perform asynchronous copy in.
+ * \param[in]     transferKind         Copy type: synchronous or asynchronous.
+ * \param[out]    timingEvent          A pointer to the H2D copy timing event to be filled in.
  *                                     If the pointer is not null, the event can further be used
  *                                     to queue a wait for this operation or to query profiling information.
  */
 template<typename ValueType>
-void copyFromDeviceBuffer(ValueType* /* hostBuffer */,
-                          DeviceBuffer<ValueType>* /* buffer */,
-                          size_t /* startingOffset */,
-                          size_t /* numValues */,
-                          const DeviceStream& /* deviceStream */,
-                          GpuApiCallBehavior /* transferKind */,
-                          CommandEvent* /* timingEvent */)
+void copyFromDeviceBuffer(ValueType* gmx_unused    hostBuffer,
+                          DeviceBuffer<ValueType>* gmx_unused buffer,
+                          size_t gmx_unused startingOffset,
+                          size_t gmx_unused   numValues,
+                          const DeviceStream& gmx_unused deviceStream,
+                          GpuApiCallBehavior gmx_unused transferKind,
+                          CommandEvent* gmx_unused timingEvent)
 {
     // SYCL-TODO
 }
@@ -146,16 +146,16 @@ void copyFromDeviceBuffer(ValueType* /* hostBuffer */,
  * Clears the device buffer asynchronously.
  *
  * \tparam        ValueType       Raw value type of the \p buffer.
- * param[in,out] buffer          Pointer to the device-side buffer
- * param[in]     startingOffset  Offset (in values) at the device-side buffer to start clearing at.
- * param[in]     numValues       Number of values to clear.
- * param[in]     deviceStream    GPU stream.
+ * \param[in,out] buffer          Pointer to the device-side buffer
+ * \param[in]     startingOffset  Offset (in values) at the device-side buffer to start clearing at.
+ * \param[in]     numValues       Number of values to clear.
+ * \param[in]     deviceStream    GPU stream.
  */
 template<typename ValueType>
-void clearDeviceBufferAsync(DeviceBuffer<ValueType>* /* buffer */,
-                            size_t /* startingOffset */,
-                            size_t /* numValues */,
-                            const DeviceStream& /* deviceStream */)
+void clearDeviceBufferAsync(DeviceBuffer<ValueType>* gmx_unused buffer,
+                            size_t gmx_unused startingOffset,
+                            size_t gmx_unused   numValues,
+                            const DeviceStream& gmx_unused deviceStream)
 {
     // SYCL-TODO
 }
@@ -164,13 +164,13 @@ void clearDeviceBufferAsync(DeviceBuffer<ValueType>* /* buffer */,
  *
  * Checks if the buffer is not nullptr and if its allocation is big enough.
  *
- * param[in] buffer        Device buffer to be checked.
- * param[in] requiredSize  Number of elements that the buffer will have to accommodate.
+ * \param[in] buffer        Device buffer to be checked.
+ * \param[in] requiredSize  Number of elements that the buffer will have to accommodate.
  *
  * \returns Whether the device buffer can be set.
  */
 template<typename T>
-static bool checkDeviceBuffer(DeviceBuffer<T> /* buffer */, int /* requiredSize */)
+static bool gmx_unused checkDeviceBuffer(DeviceBuffer<T> gmx_unused buffer, int gmx_unused requiredSize)
 {
     // SYCL-TODO
 }
@@ -186,17 +186,18 @@ using DeviceTexture = void*;
  *
  * \tparam      ValueType      Raw data type.
  *
- * param[out]  deviceBuffer   Device buffer to store data in.
- * param[in]   hostBuffer     Host buffer to get date from.
- * param[in]   numValues      Number of elements in the buffer.
- * param[in]   deviceContext  GPU device context.
+ * \param[out]  deviceBuffer   Device buffer to store data in.
+ * \param[out]  deviceTexture  New texture object
+ * \param[in]   hostBuffer     Host buffer to get date from.
+ * \param[in]   numValues      Number of elements in the buffer.
+ * \param[in]   deviceContext  GPU device context.
  */
 template<typename ValueType>
-void initParamLookupTable(DeviceBuffer<ValueType>* /* deviceBuffer */,
-                          DeviceTexture* /* deviceTexture */,
-                          const ValueType* /* hostBuffer */,
-                          int /* numValues */,
-                          const DeviceContext& /* deviceContext */)
+void initParamLookupTable(DeviceBuffer<ValueType>* gmx_unused deviceBuffer,
+                          DeviceTexture* gmx_unused deviceTexture,
+                          const ValueType* gmx_unused hostBuffer,
+                          int gmx_unused       numValues,
+                          const DeviceContext& gmx_unused deviceContext)
 {
     // SYCL-TODO
 }
@@ -205,10 +206,12 @@ void initParamLookupTable(DeviceBuffer<ValueType>* /* deviceBuffer */,
  *
  * \tparam        ValueType     Raw data type.
  *
- * param[in,out] deviceBuffer  Device buffer to store data in.
+ * \param[in,out] deviceBuffer  Device buffer to store data in.
+ * \param[in]     deviceTexture Reference to texture object
  */
 template<typename ValueType>
-void destroyParamLookupTable(DeviceBuffer<ValueType>* /* deviceBuffer */, DeviceTexture& /* deviceTexture*/)
+void destroyParamLookupTable(DeviceBuffer<ValueType>* gmx_unused deviceBuffer,
+                             DeviceTexture& gmx_unused deviceTexture)
 {
     // SYCL-TODO
 }
index f22e3208e08c29c491154122c675679c9192149f..aae05f4ad24b03c40bb898a7a815bab39a8d2c17 100644 (file)
  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
  * \inlibraryapi
  */
-#ifndef GMX_GPU_UTILS_GPUEVENTSYNCHRONIZER_OCL_H
-#define GMX_GPU_UTILS_GPUEVENTSYNCHRONIZER_OCL_H
+#ifndef GMX_GPU_UTILS_GPUEVENTSYNCHRONIZER_SYCL_H
+#define GMX_GPU_UTILS_GPUEVENTSYNCHRONIZER_SYCL_H
 
-#ifndef DOXYGEN
-
-#    include "gromacs/utility/exceptions.h"
-#    include "gromacs/utility/gmxassert.h"
+#include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/gmxassert.h"
 
+#ifndef DOXYGEN
 /*! \libinternal \brief
  * A class which allows for CPU thread to mark and wait for certain GPU stream execution point.
  * The event can be put into the stream with markEvent() and then later waited on with waitForEvent().
@@ -99,6 +98,6 @@ private:
     // SYCL-TODO
 };
 
-#endif
+#endif // !defined DOXYGEN
 
-#endif
+#endif // GMX_GPU_UTILS_GPUEVENTSYNCHRONIZER_SYCL_H