SYCL NBNXM offload support
[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 "gromacs/gpu_utils/gmxsycl.h"
48 #include "gromacs/hardware/device_management.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/stringutil.h"
51
52 #include "device_information.h"
53
54
55 bool isDeviceDetectionFunctional(std::string* errorMessage)
56 {
57     try
58     {
59         const std::vector<cl::sycl::platform> platforms = cl::sycl::platform::get_platforms();
60         // SYCL should always have the "host" platform, but just in case:
61         if (platforms.empty() && errorMessage != nullptr)
62         {
63             errorMessage->assign("No SYCL platforms found.");
64         }
65         return !platforms.empty();
66     }
67     catch (const std::exception& e)
68     {
69         if (errorMessage != nullptr)
70         {
71             errorMessage->assign(
72                     gmx::formatString("Unable to get the list of SYCL platforms: %s", e.what()));
73         }
74         return false;
75     }
76 }
77
78
79 /*!
80  * \brief Checks that device \c deviceInfo is compatible with GROMACS.
81  *
82  *  For now, only checks that the vendor is Intel and it is a GPU.
83  *
84  * \param[in]  syclDevice  The SYCL device pointer.
85  * \returns                The status enumeration value for the checked device:
86  */
87 static DeviceStatus isDeviceCompatible(const cl::sycl::device& syclDevice)
88 {
89     if (getenv("GMX_GPU_DISABLE_COMPATIBILITY_CHECK") != nullptr)
90     {
91         // Assume the device is compatible because checking has been disabled.
92         return DeviceStatus::Compatible;
93     }
94
95     if (syclDevice.get_info<cl::sycl::info::device::local_mem_type>() == cl::sycl::info::local_mem_type::none)
96     {
97         // While some kernels (leapfrog) can run without shared/local memory, this is a bad sign
98         return DeviceStatus::Incompatible;
99     }
100
101     if (syclDevice.is_host())
102     {
103         // Host device does not support subgroups or even querying for sub_group_sizes
104         return DeviceStatus::Incompatible;
105     }
106     const std::vector<size_t> supportedSubGroupSizes =
107             syclDevice.get_info<cl::sycl::info::device::sub_group_sizes>();
108     const size_t requiredSubGroupSizeForNBNXM = 8;
109     if (std::find(supportedSubGroupSizes.begin(), supportedSubGroupSizes.end(), requiredSubGroupSizeForNBNXM)
110         == supportedSubGroupSizes.end())
111     {
112         return DeviceStatus::IncompatibleClusterSize;
113     }
114
115     /* Host device can not be used, because NBNXM requires sub-groups, which are not supported.
116      * Accelerators (FPGAs and their emulators) are not supported.
117      * So, the only viable options are CPUs and GPUs. */
118     const bool forceCpu = (getenv("GMX_SYCL_FORCE_CPU") != nullptr);
119
120     if (forceCpu && syclDevice.is_cpu())
121     {
122         return DeviceStatus::Compatible;
123     }
124     else if (!forceCpu && syclDevice.is_gpu())
125     {
126         return DeviceStatus::Compatible;
127     }
128     else
129     {
130         return DeviceStatus::Incompatible;
131     }
132 }
133
134
135 /*!
136  * \brief Checks that device \c deviceInfo is sane (ie can run a kernel).
137  *
138  * Compiles and runs a dummy kernel to determine whether the given
139  * SYCL device functions properly.
140  *
141  *
142  * \param[in]  syclDevice      The device info pointer.
143  * \param[out] errorMessage    An error message related to a SYCL error.
144  * \throws     std::bad_alloc  When out of memory.
145  * \returns                    Whether the device passed sanity checks
146  */
147 static bool isDeviceFunctional(const cl::sycl::device& syclDevice, std::string* errorMessage)
148 {
149     static const int numThreads = 8;
150     try
151     {
152         cl::sycl::queue          queue(syclDevice);
153         cl::sycl::buffer<int, 1> buffer(numThreads);
154         queue.submit([&](cl::sycl::handler& cgh) {
155                  auto d_buffer = buffer.get_access<cl::sycl::access::mode::discard_write>(cgh);
156                  cl::sycl::range<1> range{ numThreads };
157                  cgh.parallel_for<class DummyKernel>(range, [=](cl::sycl::id<1> threadId) {
158                      d_buffer[threadId] = threadId.get(0);
159                  });
160              })
161                 .wait_and_throw();
162         const auto h_Buffer = buffer.get_access<cl::sycl::access::mode::read>();
163         for (int i = 0; i < numThreads; i++)
164         {
165             if (h_Buffer[i] != i)
166             {
167                 if (errorMessage != nullptr)
168                 {
169                     errorMessage->assign("Dummy kernel produced invalid values");
170                 }
171                 return false;
172             }
173         }
174     }
175     catch (const std::exception& e)
176     {
177         if (errorMessage != nullptr)
178         {
179             errorMessage->assign(
180                     gmx::formatString("Unable to run dummy kernel on device %s: %s",
181                                       syclDevice.get_info<cl::sycl::info::device::name>().c_str(),
182                                       e.what()));
183         }
184         return false;
185     }
186
187     return true;
188 }
189
190 /*!
191  * \brief Checks that device \c deviceInfo is compatible and functioning.
192  *
193  * Checks the given SYCL device for compatibility and runs a dummy kernel on it to determine
194  * whether the device functions properly.
195  *
196  *
197  * \param[in] deviceId         Device number (internal to GROMACS).
198  * \param[in] deviceInfo       The device info pointer.
199  * \returns                    The status of device.
200  */
201 static DeviceStatus checkDevice(size_t deviceId, const DeviceInformation& deviceInfo)
202 {
203
204     DeviceStatus supportStatus = isDeviceCompatible(deviceInfo.syclDevice);
205     if (supportStatus != DeviceStatus::Compatible)
206     {
207         return supportStatus;
208     }
209
210     std::string errorMessage;
211     if (!isDeviceFunctional(deviceInfo.syclDevice, &errorMessage))
212     {
213         gmx_warning("While sanity checking device #%zu, %s", deviceId, errorMessage.c_str());
214         return DeviceStatus::NonFunctional;
215     }
216
217     return DeviceStatus::Compatible;
218 }
219
220 std::vector<std::unique_ptr<DeviceInformation>> findDevices()
221 {
222     std::vector<std::unique_ptr<DeviceInformation>> deviceInfos(0);
223     const std::vector<cl::sycl::device>             devices = cl::sycl::device::get_devices();
224     deviceInfos.reserve(devices.size());
225     for (const auto& syclDevice : devices)
226     {
227         deviceInfos.emplace_back(std::make_unique<DeviceInformation>());
228
229         size_t i = deviceInfos.size() - 1;
230
231         deviceInfos[i]->id         = i;
232         deviceInfos[i]->syclDevice = syclDevice;
233         deviceInfos[i]->status     = checkDevice(i, *deviceInfos[i]);
234         deviceInfos[i]->deviceVendor =
235                 getDeviceVendor(syclDevice.get_info<cl::sycl::info::device::vendor>().c_str());
236     }
237     return deviceInfos;
238 }
239
240 void setActiveDevice(const DeviceInformation& /*deviceInfo*/) {}
241
242 void releaseDevice(DeviceInformation* /* deviceInfo */) {}
243
244 std::string getDeviceInformationString(const DeviceInformation& deviceInfo)
245 {
246     bool deviceExists = (deviceInfo.status != DeviceStatus::Nonexistent
247                          && deviceInfo.status != DeviceStatus::NonFunctional);
248
249     if (!deviceExists)
250     {
251         return gmx::formatString(
252                 "#%d: %s, status: %s", deviceInfo.id, "N/A", c_deviceStateString[deviceInfo.status]);
253     }
254     else
255     {
256         return gmx::formatString(
257                 "#%d: name: %s, vendor: %s, device version: %s, status: %s",
258                 deviceInfo.id,
259                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::name>().c_str(),
260                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::vendor>().c_str(),
261                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::version>().c_str(),
262                 c_deviceStateString[deviceInfo.status]);
263     }
264 }