93021e418f9420c975d46fe9f0b45d90d17944be
[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.is_accelerator()) // FPGAs and FPGA emulators
96     {
97         return DeviceStatus::Incompatible;
98     }
99     else
100     {
101         return DeviceStatus::Compatible;
102     }
103 }
104
105
106 /*!
107  * \brief Checks that device \c deviceInfo is sane (ie can run a kernel).
108  *
109  * Compiles and runs a dummy kernel to determine whether the given
110  * SYCL device functions properly.
111  *
112  *
113  * \param[in]  syclDevice      The device info pointer.
114  * \param[out] errorMessage    An error message related to a SYCL error.
115  * \throws     std::bad_alloc  When out of memory.
116  * \returns                    Whether the device passed sanity checks
117  */
118 static bool isDeviceFunctional(const cl::sycl::device& syclDevice, std::string* errorMessage)
119 {
120     static const int numThreads = 8;
121     try
122     {
123         cl::sycl::queue          queue(syclDevice);
124         cl::sycl::buffer<int, 1> buffer(numThreads);
125         queue.submit([&](cl::sycl::handler& cgh) {
126                  auto d_buffer = buffer.get_access<cl::sycl::access::mode::discard_write>(cgh);
127                  cl::sycl::range<1> range{ numThreads };
128                  cgh.parallel_for<class DummyKernel>(range, [=](cl::sycl::id<1> threadId) {
129                      d_buffer[threadId] = threadId.get(0);
130                  });
131              })
132                 .wait_and_throw();
133         const auto h_Buffer = buffer.get_access<cl::sycl::access::mode::read>();
134         for (int i = 0; i < numThreads; i++)
135         {
136             if (h_Buffer[i] != i)
137             {
138                 if (errorMessage != nullptr)
139                 {
140                     errorMessage->assign("Dummy kernel produced invalid values");
141                 }
142                 return false;
143             }
144         }
145     }
146     catch (const std::exception& e)
147     {
148         if (errorMessage != nullptr)
149         {
150             errorMessage->assign(
151                     gmx::formatString("Unable to run dummy kernel on device %s: %s",
152                                       syclDevice.get_info<cl::sycl::info::device::name>().c_str(),
153                                       e.what()));
154         }
155         return false;
156     }
157
158     return true;
159 }
160
161 /*!
162  * \brief Checks that device \c deviceInfo is compatible and functioning.
163  *
164  * Checks the given SYCL device for compatibility and runs a dummy kernel on it to determine
165  * whether the device functions properly.
166  *
167  *
168  * \param[in] deviceId         Device number (internal to GROMACS).
169  * \param[in] deviceInfo       The device info pointer.
170  * \returns                    The status of device.
171  */
172 static DeviceStatus checkDevice(size_t deviceId, const DeviceInformation& deviceInfo)
173 {
174
175     DeviceStatus supportStatus = isDeviceCompatible(deviceInfo.syclDevice);
176     if (supportStatus != DeviceStatus::Compatible)
177     {
178         return supportStatus;
179     }
180
181     std::string errorMessage;
182     if (!isDeviceFunctional(deviceInfo.syclDevice, &errorMessage))
183     {
184         gmx_warning("While sanity checking device #%zu, %s", deviceId, errorMessage.c_str());
185         return DeviceStatus::NonFunctional;
186     }
187
188     return DeviceStatus::Compatible;
189 }
190
191 std::vector<std::unique_ptr<DeviceInformation>> findDevices()
192 {
193     std::vector<std::unique_ptr<DeviceInformation>> deviceInfos(0);
194     const std::vector<cl::sycl::device>             devices = cl::sycl::device::get_devices();
195     deviceInfos.reserve(devices.size());
196     for (const auto& syclDevice : devices)
197     {
198         deviceInfos.emplace_back(std::make_unique<DeviceInformation>());
199
200         size_t i = deviceInfos.size() - 1;
201
202         deviceInfos[i]->id         = i;
203         deviceInfos[i]->syclDevice = syclDevice;
204         deviceInfos[i]->status     = checkDevice(i, *deviceInfos[i]);
205         deviceInfos[i]->deviceVendor =
206                 getDeviceVendor(syclDevice.get_info<cl::sycl::info::device::vendor>().c_str());
207     }
208     return deviceInfos;
209 }
210
211 void setActiveDevice(const DeviceInformation& /*deviceInfo*/) {}
212
213 void releaseDevice(DeviceInformation* /* deviceInfo */) {}
214
215 std::string getDeviceInformationString(const DeviceInformation& deviceInfo)
216 {
217     bool deviceExists = (deviceInfo.status != DeviceStatus::Nonexistent
218                          && deviceInfo.status != DeviceStatus::NonFunctional);
219
220     if (!deviceExists)
221     {
222         return gmx::formatString(
223                 "#%d: %s, status: %s", deviceInfo.id, "N/A", c_deviceStateString[deviceInfo.status]);
224     }
225     else
226     {
227         return gmx::formatString(
228                 "#%d: name: %s, vendor: %s, device version: %s, status: %s",
229                 deviceInfo.id,
230                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::name>().c_str(),
231                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::vendor>().c_str(),
232                 deviceInfo.syclDevice.get_info<cl::sycl::info::device::version>().c_str(),
233                 c_deviceStateString[deviceInfo.status]);
234     }
235 }