3b687f3660e7f382c654393aefb258737c7bc937
[alexxy/gromacs.git] / src / gromacs / hardware / device_management_sycl.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020,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  *  \brief Defines the SYCL implementations of the device management.
37  *
38  *  \author Paul Bauer <paul.bauer.q@gmail.com>
39  *  \author Erik Lindahl <erik.lindahl@gmail.com>
40  *  \author Artem Zhmurov <zhmurov@gmail.com>
41  *  \author Andrey Alekseenko <al42and@gmail.com>
42  *
43  * \ingroup module_hardware
44  */
45 #include "gmxpre.h"
46
47 #include <map>
48
49 #include "gromacs/gpu_utils/gmxsycl.h"
50 #include "gromacs/hardware/device_management.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/stringutil.h"
53
54 #include "device_information.h"
55
56
57 bool isDeviceDetectionFunctional(std::string* errorMessage)
58 {
59     try
60     {
61         const std::vector<cl::sycl::platform> platforms = cl::sycl::platform::get_platforms();
62         // SYCL should always have the "host" platform, but just in case:
63         if (platforms.empty() && errorMessage != nullptr)
64         {
65             errorMessage->assign("No SYCL platforms found.");
66         }
67         return !platforms.empty();
68     }
69     catch (const std::exception& e)
70     {
71         if (errorMessage != nullptr)
72         {
73             errorMessage->assign(
74                     gmx::formatString("Unable to get the list of SYCL platforms: %s", e.what()));
75         }
76         return false;
77     }
78 }
79
80
81 /*!
82  * \brief Checks that device \c deviceInfo is compatible with GROMACS.
83  *
84  *  For now, only checks that the vendor is Intel and it is a GPU.
85  *
86  * \param[in]  syclDevice  The SYCL device pointer.
87  * \returns                The status enumeration value for the checked device:
88  */
89 static DeviceStatus isDeviceCompatible(const cl::sycl::device& syclDevice)
90 {
91     if (getenv("GMX_GPU_DISABLE_COMPATIBILITY_CHECK") != nullptr)
92     {
93         // Assume the device is compatible because checking has been disabled.
94         return DeviceStatus::Compatible;
95     }
96
97     if (syclDevice.get_info<cl::sycl::info::device::local_mem_type>() == cl::sycl::info::local_mem_type::none)
98     {
99         // While some kernels (leapfrog) can run without shared/local memory, this is a bad sign
100         return DeviceStatus::Incompatible;
101     }
102
103     if (syclDevice.is_host())
104     {
105         // Host device does not support subgroups or even querying for sub_group_sizes
106         return DeviceStatus::Incompatible;
107     }
108
109 #if GMX_SYCL_HIPSYCL
110     /* At the time of writing:
111      * 1. SYCL NB kernels currently don't support sub_group size of 32 or 64, which are the only
112      * ones available on NVIDIA and AMD hardware, respectively. That's not a fundamental limitation,
113      * but requires porting more OpenCL code, see #3934.
114      * 2. hipSYCL does not support cl::sycl::info::device::sub_group_sizes,
115      * see https://github.com/illuhad/hipSYCL/pull/449
116      */
117     const std::vector<size_t> supportedSubGroupSizes{ warpSize };
118 #else
119     const std::vector<size_t> supportedSubGroupSizes =
120             syclDevice.get_info<cl::sycl::info::device::sub_group_sizes>();
121 #endif
122     const size_t requiredSubGroupSizeForNBNXM = 8;
123     if (std::find(supportedSubGroupSizes.begin(), supportedSubGroupSizes.end(), requiredSubGroupSizeForNBNXM)
124         == supportedSubGroupSizes.end())
125     {
126         return DeviceStatus::IncompatibleClusterSize;
127     }
128
129     /* Host device can not be used, because NBNXM requires sub-groups, which are not supported.
130      * Accelerators (FPGAs and their emulators) are not supported.
131      * So, the only viable options are CPUs and GPUs. */
132     const bool forceCpu = (getenv("GMX_SYCL_FORCE_CPU") != nullptr);
133
134     if (forceCpu && syclDevice.is_cpu())
135     {
136         return DeviceStatus::Compatible;
137     }
138     else if (!forceCpu && syclDevice.is_gpu())
139     {
140         return DeviceStatus::Compatible;
141     }
142     else
143     {
144         return DeviceStatus::Incompatible;
145     }
146 }
147
148 // Declaring the class here to avoid long unreadable name in the profiler report
149 //! \brief Class name for test kernel
150 class DummyKernel;
151
152 /*!
153  * \brief Checks that device \c deviceInfo is sane (ie can run a kernel).
154  *
155  * Compiles and runs a dummy kernel to determine whether the given
156  * SYCL device functions properly.
157  *
158  *
159  * \param[in]  syclDevice      The device info pointer.
160  * \param[out] errorMessage    An error message related to a SYCL error.
161  * \throws     std::bad_alloc  When out of memory.
162  * \returns                    Whether the device passed sanity checks
163  */
164 static bool isDeviceFunctional(const cl::sycl::device& syclDevice, std::string* errorMessage)
165 {
166     static const int numThreads = 8;
167     try
168     {
169         cl::sycl::queue          queue(syclDevice);
170         cl::sycl::buffer<int, 1> buffer(numThreads);
171         queue.submit([&](cl::sycl::handler& cgh) {
172                  auto d_buffer = buffer.get_access<cl::sycl::access::mode::discard_write>(cgh);
173                  cl::sycl::range<1> range{ numThreads };
174                  cgh.parallel_for<DummyKernel>(range, [=](cl::sycl::id<1> threadId) {
175                      d_buffer[threadId] = threadId.get(0);
176                  });
177              }).wait_and_throw();
178         const auto h_Buffer = buffer.get_access<cl::sycl::access::mode::read>();
179         for (int i = 0; i < numThreads; i++)
180         {
181             if (h_Buffer[i] != i)
182             {
183                 if (errorMessage != nullptr)
184                 {
185                     errorMessage->assign("Dummy kernel produced invalid values");
186                 }
187                 return false;
188             }
189         }
190     }
191     catch (const std::exception& e)
192     {
193         if (errorMessage != nullptr)
194         {
195             errorMessage->assign(
196                     gmx::formatString("Unable to run dummy kernel on device %s: %s",
197                                       syclDevice.get_info<cl::sycl::info::device::name>().c_str(),
198                                       e.what()));
199         }
200         return false;
201     }
202
203     return true;
204 }
205
206 /*!
207  * \brief Checks that device \c deviceInfo is compatible and functioning.
208  *
209  * Checks the given SYCL device for compatibility and runs a dummy kernel on it to determine
210  * whether the device functions properly.
211  *
212  *
213  * \param[in] deviceId         Device number (internal to GROMACS).
214  * \param[in] deviceInfo       The device info pointer.
215  * \returns                    The status of device.
216  */
217 static DeviceStatus checkDevice(size_t deviceId, const DeviceInformation& deviceInfo)
218 {
219
220     DeviceStatus supportStatus = isDeviceCompatible(deviceInfo.syclDevice);
221     if (supportStatus != DeviceStatus::Compatible)
222     {
223         return supportStatus;
224     }
225
226     std::string errorMessage;
227     if (!isDeviceFunctional(deviceInfo.syclDevice, &errorMessage))
228     {
229         gmx_warning("While sanity checking device #%zu, %s", deviceId, errorMessage.c_str());
230         return DeviceStatus::NonFunctional;
231     }
232
233     return DeviceStatus::Compatible;
234 }
235
236 /* In DPCPP, the same physical device can appear as different virtual devices provided
237  * by different backends (e.g., the same GPU can be accessible via both OpenCL and L0).
238  * Thus, using devices from two backends is more likely to be a user error than the
239  * desired behavior. In this function, we choose the backend with the most compatible
240  * devices. In case of a tie, we choose OpenCL (if present), or some arbitrary backend
241  * among those with the most devices.
242  *
243  * In hipSYCL, this problem is unlikely to manifest. It has (as of 2021-03-03) another
244  * issues: D2D copy between different backends is not allowed. We don't use D2D in
245  * SYCL yet. Additionally, hipSYCL does not implement the `sycl::platform::get_backend()`
246  * function.
247  * Thus, we only do the backend filtering with DPCPP.
248  * */
249 #if GMX_SYCL_DPCPP
250 static std::optional<cl::sycl::backend>
251 chooseBestBackend(const std::vector<std::unique_ptr<DeviceInformation>>& deviceInfos)
252 {
253     // Count the number of compatible devices per backend
254     std::map<cl::sycl::backend, int> countDevicesByBackend; // Default initialized with zeros
255     for (const auto& deviceInfo : deviceInfos)
256     {
257         if (deviceInfo->status == DeviceStatus::Compatible)
258         {
259             const cl::sycl::backend backend = deviceInfo->syclDevice.get_platform().get_backend();
260             ++countDevicesByBackend[backend];
261         }
262     }
263     // If we have devices from more than one backend...
264     if (countDevicesByBackend.size() > 1)
265     {
266         // Find backend with most devices
267         const auto backendWithMostDevices = std::max_element(
268                 countDevicesByBackend.cbegin(),
269                 countDevicesByBackend.cend(),
270                 [](const auto& kv1, const auto& kv2) { return kv1.second < kv2.second; });
271         // Count devices provided by OpenCL. Will be zero if no OpenCL devices found.
272         const int devicesInOpenCL = countDevicesByBackend[cl::sycl::backend::opencl];
273         if (devicesInOpenCL == backendWithMostDevices->second)
274         {
275             // Prefer OpenCL backend as more stable, if it has as many devices as others
276             return cl::sycl::backend::opencl;
277         }
278         else
279         {
280             // Otherwise, just return max
281             return backendWithMostDevices->first;
282         }
283     }
284     else if (countDevicesByBackend.size() == 1)
285     {
286         return countDevicesByBackend.cbegin()->first;
287     }
288     else // No devices found
289     {
290         return std::nullopt;
291     }
292 }
293 #endif
294
295 std::vector<std::unique_ptr<DeviceInformation>> findDevices()
296 {
297     std::vector<std::unique_ptr<DeviceInformation>> deviceInfos(0);
298     const std::vector<cl::sycl::device>             devices = cl::sycl::device::get_devices();
299     deviceInfos.reserve(devices.size());
300     for (const auto& syclDevice : devices)
301     {
302         deviceInfos.emplace_back(std::make_unique<DeviceInformation>());
303
304         size_t i = deviceInfos.size() - 1;
305
306         deviceInfos[i]->id         = i;
307         deviceInfos[i]->syclDevice = syclDevice;
308         deviceInfos[i]->status     = checkDevice(i, *deviceInfos[i]);
309         deviceInfos[i]->deviceVendor =
310                 getDeviceVendor(syclDevice.get_info<cl::sycl::info::device::vendor>().c_str());
311     }
312 #if GMX_SYCL_DPCPP
313     // Now, filter by the backend if we did not disable compatibility check
314     if (getenv("GMX_GPU_DISABLE_COMPATIBILITY_CHECK") == nullptr)
315     {
316         std::optional<cl::sycl::backend> preferredBackend = chooseBestBackend(deviceInfos);
317         if (preferredBackend.has_value())
318         {
319             for (auto& deviceInfo : deviceInfos)
320             {
321                 if (deviceInfo->syclDevice.get_platform().get_backend() != *preferredBackend)
322                 {
323                     deviceInfo->status = DeviceStatus::NotPreferredBackend;
324                 }
325             }
326         }
327     }
328 #endif
329     return deviceInfos;
330 }
331
332 void setActiveDevice(const DeviceInformation& /*deviceInfo*/) {}
333
334 void releaseDevice(DeviceInformation* /* deviceInfo */) {}
335
336 std::string getDeviceInformationString(const DeviceInformation& deviceInfo)
337 {
338     bool deviceExists = (deviceInfo.status != DeviceStatus::Nonexistent
339                          && deviceInfo.status != DeviceStatus::NonFunctional);
340
341     if (!deviceExists)
342     {
343         return gmx::formatString(
344                 "#%d: %s, status: %s", deviceInfo.id, "N/A", c_deviceStateString[deviceInfo.status]);
345     }
346     else
347     {
348         return gmx::formatString(
349                 "#%d: name: %s, vendor: %s, device version: %s, status: %s",
350                 deviceInfo.id,
351                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::name>().c_str(),
352                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::vendor>().c_str(),
353                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::version>().c_str(),
354                 c_deviceStateString[deviceInfo.status]);
355     }
356 }