Draft: Fix compilation
authorAndrey Alekseenko <al42and@gmail.com>
Tue, 22 Sep 2020 12:50:41 +0000 (12:50 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 22 Sep 2020 12:50:41 +0000 (12:50 +0000)
Also, make ctor etc always defined for DeviceInformation. Otherwise it's
just ugly.

17 files changed:
cmake/gmxManageSYCL.cmake
docs/user-guide/environment-variables.rst
src/gromacs/CMakeLists.txt
src/gromacs/gpu_utils/gmxsycl.h [new file with mode: 0644]
src/gromacs/gpu_utils/tests/CMakeLists.txt
src/gromacs/gpu_utils/tests/device_stream_manager.cpp
src/gromacs/hardware/CMakeLists.txt
src/gromacs/hardware/detecthardware.cpp
src/gromacs/hardware/device_information.h
src/gromacs/hardware/device_management.cpp
src/gromacs/hardware/device_management_common.cpp
src/gromacs/hardware/device_management_ocl.cpp
src/gromacs/hardware/device_management_sycl.cpp [new file with mode: 0644]
src/gromacs/hardware/tests/CMakeLists.txt
src/gromacs/hardware/tests/device_management.cpp
src/gromacs/taskassignment/CMakeLists.txt
src/testutils/TestMacros.cmake

index e19cafb6842f41d2111aa515fa293ffb6b987b77..57c1338d5cd7a73e57ad21e8fd3bee728b7b53be 100644 (file)
@@ -65,7 +65,7 @@ if(CMAKE_CXX_COMPILER MATCHES "dpcpp")
     gmx_check_source_compiles_with_flags(
         "int main() { return 0; }"
         "-fno-sycl"
-       "CXX"
+        "CXX"
         DISABLE_SYCL_CXX_FLAGS_RESULT)
 
     if(DISABLE_SYCL_CXX_FLAGS_RESULT)
index 9868a443a21e4d850ce6d0cf4577e93a3ba46c27..8d288eed43bda79949bcc8e2ed8e42662884a4f7 100644 (file)
@@ -253,6 +253,10 @@ Performance and Run Control
         runtime permits this variable to be different for different ranks. Cannot be used
         in conjunction with ``mdrun -gputasks``. Has all the same requirements as ``mdrun -gputasks``.
 
+``GMX_GPU_DISABLE_COMPATIBILITY_CHECK``
+        Disables the hardware compatibility check in OpenCL and SYCL. Useful for developers
+        and allows testing the OpenCL/SYCL kernels on non-supported platforms without source code modification.
+
 ``GMX_IGNORE_FSYNC_FAILURE_ENV``
         allow :ref:`gmx mdrun` to continue even if
         a file is missing.
@@ -475,11 +479,6 @@ compilation of OpenCL kernels, but they are also used in device selection.
         override |Gromacs| default behavior, or if you want to test
         your own kernels.
 
-``GMX_OCL_DISABLE_COMPATIBILITY_CHECK``
-        Disables the hardware compatibility check. Useful for developers
-        and allows testing the OpenCL kernels on non-supported platforms
-        (like Intel iGPUs) without source code modification.
-
 ``GMX_OCL_SHOW_DIAGNOSTICS``
         Use Intel OpenCL extension to show additional runtime performance
         diagnostics.
index ec0422488435cd7115fa6f65af02e70abd36b985..fc791d70334042c26e6459524e5ca3e75f372de0 100644 (file)
@@ -299,7 +299,7 @@ target_link_libraries(libgromacs
                       ${GMX_COMMON_LIBRARIES}
                       ${FFT_LIBRARIES} ${LINEAR_ALGEBRA_LIBRARIES}
                       ${THREAD_LIB} ${GMX_SHARED_LINKER_FLAGS}
-                     ${SYCL_CXX_FLAGS}
+                      ${SYCL_CXX_FLAGS}
                       ${OpenCL_LIBRARIES}
                       $<$<PLATFORM_ID:SunOS>:socket>
                       PUBLIC
diff --git a/src/gromacs/gpu_utils/gmxsycl.h b/src/gromacs/gpu_utils/gmxsycl.h
new file mode 100644 (file)
index 0000000..6c1b635
--- /dev/null
@@ -0,0 +1,61 @@
+/*
+ * 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.
+ */
+/*! \libinternal \file
+ * \brief
+ * Wraps the complexity of including SYCL in GROMACS.
+ *
+ * SYCL headers use symbol DIM as a template parameter, which gets broken by macro DIM defined
+ * in gromacs/math/vectypes.h. Here, we include the SYCL header while temporary undefining this macro.
+ *
+ * \inlibraryapi
+ */
+
+#ifndef GMX_GPU_UTILS_GMXSYCL_H
+#define GMX_GPU_UTILS_GMXSYCL_H
+
+
+#ifdef DIM
+#    if DIM != 3
+#        error "The workaround here assumes we use DIM=3."
+#    else
+#        undef DIM
+#        include <CL/sycl.hpp>
+#        define DIM 3
+#    endif
+#else
+#    include <CL/sycl.hpp>
+#endif
+
+#endif
index 41407185195bc021d0846b4d3c813bdaa395991b..ed5de307d151d490e58ba405a039be7b8781e9fd 100644 (file)
@@ -57,4 +57,7 @@ gmx_add_unit_test(GpuUtilsUnitTests gpu_utils-test HARDWARE_DETECTION
 
     NON_GPU_CPP_SOURCE_FILES
         devicetransfers.cpp
+
+    SYCL_CPP_SOURCE_FILES  # SYCL-TODO: proper test
+        devicetransfers.cpp
     )
index d01aff20ead9771860aa31cdf0f19893cafe2eea..06c03fe8b1f72a30eb9b3694d9dfaf4cd0dbbeef 100644 (file)
@@ -73,7 +73,7 @@ const EnumerationArray<DeviceStreamType, std::string> c_deviceStreamNames = {
 
 /*! \brief Non-GPU builds return nullptr instead of streams,
  * so we have to expect that in such build configurations. */
-constexpr bool c_canExpectValidStreams = (GMX_GPU != 0);
+constexpr bool c_canExpectValidStreams = (GMX_GPU != 0 && !GMX_GPU_SYCL); // SYCL-TODO
 
 //! Helper function to implement readable testing
 void expectValidStreams(DeviceStreamManager* manager, std::initializer_list<DeviceStreamType> types)
index c93767f790ee008717e341c6b0a9bb0023e641b7..d51921d222adf80679d3fb1a1a17d4e5f98f5315 100644 (file)
@@ -49,6 +49,16 @@ elseif(GMX_GPU_CUDA)
     gmx_add_libgromacs_sources(
         device_management.cu
         )
+elseif(GMX_GPU_SYCL)
+    gmx_add_libgromacs_sources(
+      device_management_sycl.cpp
+    )
+    _gmx_add_files_to_property(SYCL_SOURCES
+      device_management_sycl.cpp
+      # Must add these files so they can include device_information.h
+      device_management.cpp
+      device_management_common.cpp
+      detecthardware.cpp)
 else()
     gmx_add_libgromacs_sources(
         device_management.cpp
index b1bacbeebef17ec1f17ac28d37d811c0f4123aff..3bcdac4835f2d52b5db9855e7aaab9109b213f40 100644 (file)
@@ -131,7 +131,7 @@ static void gmx_detect_gpus(const gmx::MDLogger&             mdlog,
     /* The OpenCL support requires us to run detection on all ranks.
      * With CUDA we don't need to, and prefer to detect on one rank
      * and send the information to the other ranks over MPI. */
-    constexpr bool allRanksMustDetectGpus = (GMX_GPU_OPENCL != 0);
+    constexpr bool allRanksMustDetectGpus = (GMX_GPU_OPENCL != 0 || GMX_GPU_SYCL != 0);
     bool           gpusCanBeDetected      = false;
     if (isMasterRankOfPhysicalNode || allRanksMustDetectGpus)
     {
index e3005e1d98ae2629ff75a534e8cad1f7af31cafc..fd1a806b886ce145263ea72f8a33b6f30b2acbe6 100644 (file)
 #if GMX_GPU_OPENCL
 #    include "gromacs/gpu_utils/gmxopencl.h"
 #endif
+
+#if GMX_GPU_SYCL
+#    include "gromacs/gpu_utils/gmxsycl.h"
+#endif
+
 #include "gromacs/utility/enumerationhelpers.h"
 
 //! Constant used to help minimize preprocessed code
 static constexpr bool c_binarySupportsGpus = (GMX_GPU != 0);
+static constexpr bool c_canSerializeDeviceInformation =
+        (!GMX_GPU_OPENCL && !GMX_GPU_SYCL); /*NOLINT(misc-redundant-expression)*/
 
 //! Possible results of the GPU detection/check.
 enum class DeviceStatus : int
@@ -147,6 +154,8 @@ struct DeviceInformation
     DeviceVendor   deviceVendor;        //!< Device vendor.
     size_t         maxWorkItemSizes[3]; //!< Workgroup size limits (CL_DEVICE_MAX_WORK_ITEM_SIZES).
     size_t         maxWorkGroupSize;    //!< Workgroup total size limit (CL_DEVICE_MAX_WORK_GROUP_SIZE).
+#elif GMX_GPU_SYCL
+    cl::sycl::device syclDevice;
 #endif
 };
 
index 873d1258c153ff57ccba5f06443927662a8e3b44..87384c78764d6eafa0ab46f9a7fefc6e7005b04e 100644 (file)
@@ -49,7 +49,6 @@
 
 #include "device_management.h"
 
-#include "gromacs/gpu_utils/gputraits.h"
 #include "gromacs/utility/fatalerror.h"
 
 #include "device_information.h"
index 3d043234321f3df055a78cd76fb54da9ac709711..5b2e5af5d328db02eae9cb01865eca2660a0561d 100644 (file)
@@ -108,6 +108,8 @@ std::string getDeviceCompatibilityDescription(const std::vector<std::unique_ptr<
 void serializeDeviceInformations(const std::vector<std::unique_ptr<DeviceInformation>>& deviceInfoList,
                                  gmx::ISerializer*                                      serializer)
 {
+    GMX_RELEASE_ASSERT(c_canSerializeDeviceInformation,
+                       "DeviceInformation for OpenCL/SYCL can not be serialized");
     int numDevices = deviceInfoList.size();
     serializer->doInt(&numDevices);
     for (auto& deviceInfo : deviceInfoList)
@@ -118,6 +120,8 @@ void serializeDeviceInformations(const std::vector<std::unique_ptr<DeviceInforma
 
 std::vector<std::unique_ptr<DeviceInformation>> deserializeDeviceInformations(gmx::ISerializer* serializer)
 {
+    GMX_RELEASE_ASSERT(c_canSerializeDeviceInformation,
+                       "DeviceInformation for OpenCL/SYCL can not be deserialized");
     int numDevices = 0;
     serializer->doInt(&numDevices);
     std::vector<std::unique_ptr<DeviceInformation>> deviceInfoList(numDevices);
index 5380383fcfd68a1e2778bb0d838ee64fd2572151..14d211089996c008c774cbc5d4582cfa4222d6b3 100644 (file)
@@ -152,11 +152,19 @@ static bool runningOnCompatibleHWForNvidia(const DeviceInformation& deviceInfo)
  */
 static DeviceStatus isDeviceFunctional(const DeviceInformation& deviceInfo)
 {
-    if (getenv("GMX_OCL_DISABLE_COMPATIBILITY_CHECK") != nullptr)
+    if (getenv("GMX_GPU_DISABLE_COMPATIBILITY_CHECK") != nullptr)
     {
         // Assume the device is compatible because checking has been disabled.
         return DeviceStatus::Compatible;
     }
+    if (getenv("GMX_OCL_DISABLE_COMPATIBILITY_CHECK") != nullptr)
+    {
+        fprintf(stderr,
+                "Environment variable GMX_OCL_DISABLE_COMPATIBILITY_CHECK is deprecated and will "
+                "be removed in release 2022. Please use GMX_GPU_DISABLE_COMPATIBILITY_CHECK "
+                "instead.\n");
+        return DeviceStatus::Compatible;
+    }
 
     // OpenCL device version check, ensure >= REQUIRED_OPENCL_MIN_VERSION
     constexpr unsigned int minVersionMajor = REQUIRED_OPENCL_MIN_VERSION_MAJOR;
diff --git a/src/gromacs/hardware/device_management_sycl.cpp b/src/gromacs/hardware/device_management_sycl.cpp
new file mode 100644 (file)
index 0000000..8a7e74e
--- /dev/null
@@ -0,0 +1,231 @@
+/*
+ * 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 Defines the SYCL implementations of the device management.
+ *
+ *  \author Paul Bauer <paul.bauer.q@gmail.com>
+ *  \author Erik Lindahl <erik.lindahl@gmail.com>
+ *  \author Artem Zhmurov <zhmurov@gmail.com>
+ *  \author Andrey Alekseenko <al42and@gmail.com>
+ *
+ * \ingroup module_hardware
+ */
+#include "gmxpre.h"
+
+#include "gromacs/gpu_utils/gmxsycl.h"
+#include "gromacs/hardware/device_management.h"
+#include "gromacs/utility/fatalerror.h"
+#include "gromacs/utility/stringutil.h"
+
+#include "device_information.h"
+
+
+bool isDeviceDetectionFunctional(std::string* errorMessage)
+{
+    try
+    {
+        const std::vector<cl::sycl::platform> platforms = cl::sycl::platform::get_platforms();
+        // SYCL should always have the "host" platform, but just in case:
+        if (platforms.empty() && errorMessage != nullptr)
+        {
+            errorMessage->assign("No SYCL platforms found.");
+        }
+        return !platforms.empty();
+    }
+    catch (const std::exception& e)
+    {
+        if (errorMessage != nullptr)
+        {
+            errorMessage->assign(
+                    gmx::formatString("Unable to get the list of SYCL platforms: %s", e.what()));
+        }
+        return false;
+    }
+}
+
+
+/*!
+ * \brief Checks that device \c deviceInfo is compatible with GROMACS.
+ *
+ *  For now, only checks that the vendor is Intel and it is a GPU.
+ *
+ * \param[in]  syclDevice  The SYCL device pointer.
+ * \returns                The status enumeration value for the checked device:
+ */
+static DeviceStatus isDeviceCompatible(const cl::sycl::device& syclDevice)
+{
+    if (getenv("GMX_GPU_DISABLE_COMPATIBILITY_CHECK") != nullptr)
+    {
+        // Assume the device is compatible because checking has been disabled.
+        return DeviceStatus::Compatible;
+    }
+
+    if (syclDevice.is_accelerator()) // FPGAs and FPGA emulators
+    {
+        return DeviceStatus::Incompatible;
+    }
+    else
+    {
+        return DeviceStatus::Compatible;
+    }
+}
+
+
+/*!
+ * \brief Checks that device \c deviceInfo is sane (ie can run a kernel).
+ *
+ * Compiles and runs a dummy kernel to determine whether the given
+ * SYCL device functions properly.
+ *
+ *
+ * \param[in]  syclDevice      The device info pointer.
+ * \param[out] errorMessage    An error message related to a SYCL error.
+ * \throws     std::bad_alloc  When out of memory.
+ * \returns                    Whether the device passed sanity checks
+ */
+static bool isDeviceFunctional(const cl::sycl::device& syclDevice, std::string* errorMessage)
+{
+    static const int numThreads = 8;
+    cl::sycl::queue  queue;
+    try
+    {
+        queue = cl::sycl::queue(syclDevice);
+        cl::sycl::buffer<int, 1> buffer(numThreads);
+        queue.submit([&](cl::sycl::handler& cgh) {
+                 auto d_buffer = buffer.get_access<cl::sycl::access::mode::discard_write>(cgh);
+                 cgh.parallel_for<class DummyKernel>(numThreads, [=](cl::sycl::id<1> threadId) {
+                     d_buffer[threadId] = threadId.get(0);
+                 });
+             })
+                .wait_and_throw();
+        const auto h_Buffer = buffer.get_access<cl::sycl::access::mode::read>();
+        for (int i = 0; i < numThreads; i++)
+        {
+            if (h_Buffer[i] != i)
+            {
+                if (errorMessage != nullptr)
+                {
+                    errorMessage->assign("Dummy kernel produced invalid values");
+                }
+                return false;
+            }
+        }
+    }
+    catch (const std::exception& e)
+    {
+        if (errorMessage != nullptr)
+        {
+            errorMessage->assign(gmx::formatString(
+                    "Unable to run dummy kernel on device %s: %s",
+                    syclDevice.get_info<cl::sycl::info::device::name>().c_str(), e.what()));
+        }
+        return false;
+    }
+
+    return true;
+}
+
+/*!
+ * \brief Checks that device \c deviceInfo is compatible and functioning.
+ *
+ * Checks the given SYCL device for compatibility and runs a dummy kernel on it to determine
+ * whether the device functions properly.
+ *
+ *
+ * \param[in] deviceId         Device number (internal to GROMACS).
+ * \param[in] deviceInfo       The device info pointer.
+ * \returns                    The status of device.
+ */
+static DeviceStatus checkDevice(size_t deviceId, const DeviceInformation& deviceInfo)
+{
+
+    DeviceStatus supportStatus = isDeviceCompatible(deviceInfo.syclDevice);
+    if (supportStatus != DeviceStatus::Compatible)
+    {
+        return supportStatus;
+    }
+
+    std::string errorMessage;
+    if (!isDeviceFunctional(deviceInfo.syclDevice, &errorMessage))
+    {
+        gmx_warning("While sanity checking device #%zu, %s", deviceId, errorMessage.c_str());
+        return DeviceStatus::NonFunctional;
+    }
+
+    return DeviceStatus::Compatible;
+}
+
+std::vector<std::unique_ptr<DeviceInformation>> findDevices()
+{
+    std::vector<std::unique_ptr<DeviceInformation>> deviceInfos(0);
+    const std::vector<cl::sycl::device>             devices = cl::sycl::device::get_devices();
+    deviceInfos.reserve(devices.size());
+    for (const auto& syclDevice : devices)
+    {
+        deviceInfos.emplace_back(std::make_unique<DeviceInformation>());
+
+        size_t i = deviceInfos.size() - 1;
+
+        deviceInfos[i]->id         = i;
+        deviceInfos[i]->syclDevice = syclDevice;
+        deviceInfos[i]->status     = checkDevice(i, *deviceInfos[i]);
+    }
+    return deviceInfos;
+}
+
+void setActiveDevice(const DeviceInformation& /*deviceInfo*/) {}
+
+void releaseDevice(DeviceInformation* /* deviceInfo */) {}
+
+std::string getDeviceInformationString(const DeviceInformation& deviceInfo)
+{
+    bool deviceExists = (deviceInfo.status != DeviceStatus::Nonexistent
+                         && deviceInfo.status != DeviceStatus::NonFunctional);
+
+    if (!deviceExists)
+    {
+        return gmx::formatString("#%d: %s, status: %s", deviceInfo.id, "N/A",
+                                 c_deviceStateString[deviceInfo.status]);
+    }
+    else
+    {
+        return gmx::formatString(
+                "#%d: name: %s, vendor: %s, device version: %s, status: %s", deviceInfo.id,
+                deviceInfo.syclDevice.get_info<cl::sycl::info::device::name>().c_str(),
+                deviceInfo.syclDevice.get_info<cl::sycl::info::device::vendor>().c_str(),
+                deviceInfo.syclDevice.get_info<cl::sycl::info::device::version>().c_str(),
+                c_deviceStateString[deviceInfo.status]);
+    }
+}
index 8479ac5399354265aeadfe9074b3d77e0a0fcc1b..5bd154009ffd2d3ed193cfe56f7b91e083f762d3 100644 (file)
@@ -35,6 +35,7 @@
 gmx_add_unit_test(HardwareUnitTests hardware-test
     CPP_SOURCE_FILES
         cpuinfo.cpp
-        device_management.cpp
         hardwaretopology.cpp
+    GPU_CPP_SOURCE_FILES
+        device_management.cpp
         )
index cfbcd134fbe545b18fe9763a277da7a3479de7c2..88d61cffca1462da6f4f790ae247b970b82a3fb1 100644 (file)
@@ -58,7 +58,7 @@ namespace
 
 TEST(DevicesManagerTest, Serialization)
 {
-    if (canPerformDeviceDetection(nullptr))
+    if (canPerformDeviceDetection(nullptr) && c_canSerializeDeviceInformation)
     {
         std::vector<std::unique_ptr<DeviceInformation>> deviceInfoListIn = findDevices();
         gmx::InMemorySerializer                         writer;
index 9afd05e13b5e227d1f0f0c372510ff22f7c56cbd..1410e3188f997f5f9add86e081dc4ca56868c888 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2017,2019, by the GROMACS development team, led by
+# Copyright (c) 2017,2019,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.
@@ -42,6 +42,10 @@ gmx_add_libgromacs_sources(
     usergpuids.cpp
     )
 
+if (GMX_GPU_SYCL)
+    _gmx_add_files_to_property(SYCL_SOURCES usergpuids.cpp)
+endif()
+
 if (BUILD_TESTING)
     add_subdirectory(tests)
 endif()
index a9fb476a567abcb23ddfc4419308a66f82190f29..6e7586654a699f11590ab9be7532ca8ce45e7792 100644 (file)
@@ -73,8 +73,10 @@ endfunction ()
 #     All the normal CUDA .cu source files
 #   OPENCL_CPP_SOURCE_FILES   file1.cpp file2.cpp ...
 #     All the other C++ .cpp source files needed only with OpenCL
+#   SYCL_CPP_SOURCE_FILES   file1.cpp file2.cpp ...
+#     All the C++ .cpp source files needed only with SYCL
 #   NON_GPU_CPP_SOURCE_FILES  file1.cpp file2.cpp ...
-#     All the other C++ .cpp source files needed only with neither OpenCL nor CUDA
+#     All the other C++ .cpp source files needed only with neither OpenCL nor CUDA nor SYCL
 function (gmx_add_gtest_executable EXENAME)
     if (GMX_BUILD_UNITTESTS AND BUILD_TESTING)
         set(_options MPI HARDWARE_DETECTION)
@@ -83,6 +85,7 @@ function (gmx_add_gtest_executable EXENAME)
             CUDA_CU_SOURCE_FILES
             GPU_CPP_SOURCE_FILES
             OPENCL_CPP_SOURCE_FILES
+            SYCL_CPP_SOURCE_FILES
             NON_GPU_CPP_SOURCE_FILES
             )
         cmake_parse_arguments(ARG "${_options}" "" "${_multi_value_keywords}" ${ARGN})
@@ -140,6 +143,13 @@ function (gmx_add_gtest_executable EXENAME)
             if(ARG_OPENCL_CPP_SOURCE_FILES OR ARG_GPU_CPP_SOURCE_FILES)
                 target_link_libraries(${EXENAME} PRIVATE ${OpenCL_LIBRARIES})
             endif()
+        elseif (GMX_GPU_SYCL)
+            target_sources(${EXENAME} PRIVATE ${ARG_SYCL_CPP_SOURCE_FILES} ${ARG_GPU_CPP_SOURCE_FILES})
+            set_source_files_properties(${ARG_GPU_CPP_SOURCE_FILES} PROPERTIES COMPILE_FLAGS "${SYCL_CXX_FLAGS}")
+            set_source_files_properties(${ARG_SYCL_CPP_SOURCE_FILES} PROPERTIES COMPILE_FLAGS "${SYCL_CXX_FLAGS}")
+            if(ARG_SYCL_CPP_SOURCE_FILES OR ARG_GPU_CPP_SOURCE_FILES)
+                target_link_libraries(${EXENAME} PRIVATE ${SYCL_CXX_FLAGS})
+            endif()
         else()
             target_sources(${EXENAME} PRIVATE ${ARG_NON_GPU_CPP_SOURCE_FILES} ${ARG_GPU_CPP_SOURCE_FILES})
         endif()