b59b51e59a10af1a1226ed879123207bed8ac310
[alexxy/gromacs.git] / src / gromacs / hardware / device_management.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016, by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  *  \brief Defines the CUDA implementations of the device management.
38  *
39  *  \author Anca Hamuraru <anca@streamcomputing.eu>
40  *  \author Dimitrios Karkoulis <dimitris.karkoulis@gmail.com>
41  *  \author Teemu Virolainen <teemu@streamcomputing.eu>
42  *  \author Mark Abraham <mark.j.abraham@gmail.com>
43  *  \author Szilárd Páll <pall.szilard@gmail.com>
44  *  \author Artem Zhmurov <zhmurov@gmail.com>
45  *
46  * \ingroup module_hardware
47  */
48 #include "gmxpre.h"
49
50 #include "device_management.h"
51
52 #include <assert.h>
53
54 #include "gromacs/gpu_utils/cudautils.cuh"
55 #include "gromacs/gpu_utils/device_context.h"
56 #include "gromacs/gpu_utils/device_stream.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/programcontext.h"
59 #include "gromacs/utility/smalloc.h"
60 #include "gromacs/utility/stringutil.h"
61
62 #include "device_information.h"
63
64 /*! \internal \brief
65  * Max number of devices supported by CUDA (for consistency checking).
66  *
67  * In reality it is 16 with CUDA <=v5.0, but let's stay on the safe side.
68  */
69 static int c_cudaMaxDeviceCount = 32;
70
71 /** Dummy kernel used for sanity checking. */
72 static __global__ void dummy_kernel(void) {}
73
74 static cudaError_t checkCompiledTargetCompatibility(int deviceId, const cudaDeviceProp& deviceProp)
75 {
76     cudaFuncAttributes attributes;
77     cudaError_t        stat = cudaFuncGetAttributes(&attributes, dummy_kernel);
78
79     if (cudaErrorInvalidDeviceFunction == stat)
80     {
81         fprintf(stderr,
82                 "\nWARNING: The %s binary does not include support for the CUDA architecture of "
83                 "the GPU ID #%d (compute capability %d.%d) detected during detection. "
84                 "By default, GROMACS supports all architectures of compute "
85                 "capability >= 3.0, so your GPU "
86                 "might be rare, or some architectures were disabled in the build. \n"
87                 "Consult the install guide for how to use the GMX_CUDA_TARGET_SM and "
88                 "GMX_CUDA_TARGET_COMPUTE CMake variables to add this architecture. \n",
89                 gmx::getProgramContext().displayName(), deviceId, deviceProp.major, deviceProp.minor);
90     }
91
92     return stat;
93 }
94
95 /*!
96  * \brief Runs GPU sanity checks.
97  *
98  * Runs a series of checks to determine that the given GPU and underlying CUDA
99  * driver/runtime functions properly.
100  *
101  * \todo Currently we do not make a distinction between the type of errors
102  *       that can appear during functionality checks. This needs to be improved,
103  *       e.g if the dummy test kernel fails to execute with a "device busy message"
104  *       we should appropriately report that the device is busy instead of NonFunctional.
105  *
106  * \todo Introduce errors codes and handle errors more smoothly.
107  *
108  *
109  * \param[in]  deviceInfo  Device information on the device to check.
110  * \returns                The status enumeration value for the checked device:
111  */
112 static DeviceStatus isDeviceFunctional(const DeviceInformation& deviceInfo)
113 {
114     cudaError_t cu_err;
115
116     /* both major & minor is 9999 if no CUDA capable devices are present */
117     if (deviceInfo.prop.major == 9999 && deviceInfo.prop.minor == 9999)
118     {
119         return DeviceStatus::NonFunctional;
120     }
121     /* we don't care about emulation mode */
122     if (deviceInfo.prop.major == 0)
123     {
124         return DeviceStatus::NonFunctional;
125     }
126
127     cu_err = cudaSetDevice(deviceInfo.id);
128     if (cu_err != cudaSuccess)
129     {
130         fprintf(stderr, "Error %d while switching to device #%d: %s\n", cu_err, deviceInfo.id,
131                 cudaGetErrorString(cu_err));
132         return DeviceStatus::NonFunctional;
133     }
134
135     cu_err = checkCompiledTargetCompatibility(deviceInfo.id, deviceInfo.prop);
136     // Avoid triggering an error if GPU devices are in exclusive or prohibited mode;
137     // it is enough to check for cudaErrorDevicesUnavailable only here because
138     // if we encounter it that will happen in cudaFuncGetAttributes in the above function.
139     if (cu_err == cudaErrorDevicesUnavailable)
140     {
141         return DeviceStatus::Unavailable;
142     }
143     else if (cu_err != cudaSuccess)
144     {
145         return DeviceStatus::NonFunctional;
146     }
147
148     /* try to execute a dummy kernel */
149     try
150     {
151         KernelLaunchConfig config;
152         config.blockSize[0]                = 512;
153         const auto          dummyArguments = prepareGpuKernelArguments(dummy_kernel, config);
154         const DeviceContext deviceContext(deviceInfo);
155         const DeviceStream  deviceStream(deviceContext, DeviceStreamPriority::Normal, false);
156         launchGpuKernel(dummy_kernel, config, deviceStream, nullptr, "Dummy kernel", dummyArguments);
157     }
158     catch (gmx::GromacsException& ex)
159     {
160         // launchGpuKernel error is not fatal and should continue with marking the device bad
161         fprintf(stderr,
162                 "Error occurred while running dummy kernel sanity check on device #%d:\n %s\n",
163                 deviceInfo.id, formatExceptionMessageToString(ex).c_str());
164         return DeviceStatus::NonFunctional;
165     }
166
167     if (cudaDeviceSynchronize() != cudaSuccess)
168     {
169         return DeviceStatus::NonFunctional;
170     }
171
172     cu_err = cudaDeviceReset();
173     CU_RET_ERR(cu_err, "cudaDeviceReset failed");
174
175     return DeviceStatus::Compatible;
176 }
177
178 /*! \brief Returns true if the gpu characterized by the device properties is supported
179  *         by the native gpu acceleration.
180  *
181  * \param[in] deviceProperties  The CUDA device properties of the gpus to test.
182  * \returns                     True if the GPU properties passed indicate a compatible
183  *                              GPU, otherwise false.
184  */
185 static bool isDeviceGenerationSupported(const cudaDeviceProp& deviceProperties)
186 {
187     return (deviceProperties.major >= 3);
188 }
189
190 /*! \brief Checks if a GPU with a given ID is supported by the native GROMACS acceleration.
191  *
192  *  Returns a status value which indicates compatibility or one of the following
193  *  errors: incompatibility or insanity (=unexpected behavior).
194  *
195  *  As the error handling only permits returning the state of the GPU, this function
196  *  does not clear the CUDA runtime API status allowing the caller to inspect the error
197  *  upon return. Note that this also means it is the caller's responsibility to
198  *  reset the CUDA runtime state.
199  *
200  *  \param[in]  deviceInfo The device information on the device to check.
201  *  \returns               the status of the requested device
202  */
203 static DeviceStatus checkDeviceStatus(const DeviceInformation& deviceInfo)
204 {
205     if (!isDeviceGenerationSupported(deviceInfo.prop))
206     {
207         return DeviceStatus::Incompatible;
208     }
209     return isDeviceFunctional(deviceInfo);
210 }
211
212 bool isDeviceDetectionFunctional(std::string* errorMessage)
213 {
214     cudaError_t stat;
215     int         driverVersion = -1;
216     stat                      = cudaDriverGetVersion(&driverVersion);
217     GMX_ASSERT(stat != cudaErrorInvalidValue,
218                "An impossible null pointer was passed to cudaDriverGetVersion");
219     GMX_RELEASE_ASSERT(
220             stat == cudaSuccess,
221             gmx::formatString("An unexpected value was returned from cudaDriverGetVersion %s: %s",
222                               cudaGetErrorName(stat), cudaGetErrorString(stat))
223                     .c_str());
224     bool foundDriver = (driverVersion > 0);
225     if (!foundDriver)
226     {
227         // Can't detect GPUs if there is no driver
228         if (errorMessage != nullptr)
229         {
230             errorMessage->assign("No valid CUDA driver found");
231         }
232         return false;
233     }
234
235     int numDevices;
236     stat = cudaGetDeviceCount(&numDevices);
237     if (stat != cudaSuccess)
238     {
239         if (errorMessage != nullptr)
240         {
241             /* cudaGetDeviceCount failed which means that there is
242              * something wrong with the machine: driver-runtime
243              * mismatch, all GPUs being busy in exclusive mode,
244              * invalid CUDA_VISIBLE_DEVICES, or some other condition
245              * which should result in GROMACS issuing at least a
246              * warning. */
247             errorMessage->assign(cudaGetErrorString(stat));
248         }
249
250         // Consume the error now that we have prepared to handle
251         // it. This stops it reappearing next time we check for
252         // errors. Note that if CUDA_VISIBLE_DEVICES does not contain
253         // valid devices, then cudaGetLastError returns the
254         // (undocumented) cudaErrorNoDevice, but this should not be a
255         // problem as there should be no future CUDA API calls.
256         // NVIDIA bug report #2038718 has been filed.
257         cudaGetLastError();
258         // Can't detect GPUs
259         return false;
260     }
261
262     // We don't actually use numDevices here, that's not the job of
263     // this function.
264     return true;
265 }
266
267 std::vector<std::unique_ptr<DeviceInformation>> findDevices()
268 {
269     int         numDevices;
270     cudaError_t stat = cudaGetDeviceCount(&numDevices);
271     if (stat != cudaSuccess)
272     {
273         GMX_THROW(gmx::InternalError(
274                 "Invalid call of findDevices() when CUDA API returned an error, perhaps "
275                 "canPerformDeviceDetection() was not called appropriately beforehand."));
276     }
277
278     /* things might go horribly wrong if cudart is not compatible with the driver */
279     numDevices = std::min(numDevices, c_cudaMaxDeviceCount);
280
281     // We expect to start device support/sanity checks with a clean runtime error state
282     gmx::ensureNoPendingCudaError("");
283
284     std::vector<std::unique_ptr<DeviceInformation>> deviceInfoList(numDevices);
285     for (int i = 0; i < numDevices; i++)
286     {
287         cudaDeviceProp prop;
288         memset(&prop, 0, sizeof(cudaDeviceProp));
289         stat = cudaGetDeviceProperties(&prop, i);
290
291         deviceInfoList[i]       = std::make_unique<DeviceInformation>();
292         deviceInfoList[i]->id   = i;
293         deviceInfoList[i]->prop = prop;
294
295         const DeviceStatus checkResult = (stat != cudaSuccess) ? DeviceStatus::NonFunctional
296                                                                : checkDeviceStatus(*deviceInfoList[i]);
297
298         deviceInfoList[i]->status = checkResult;
299
300         if (checkResult != DeviceStatus::Compatible)
301         {
302             // TODO:
303             //  - we inspect the CUDA API state to retrieve and record any
304             //    errors that occurred during is_gmx_supported_gpu_id() here,
305             //    but this would be more elegant done within is_gmx_supported_gpu_id()
306             //    and only return a string with the error if one was encountered.
307             //  - we'll be reporting without rank information which is not ideal.
308             //  - we'll end up warning also in cases where users would already
309             //    get an error before mdrun aborts.
310             //
311             // Here we also clear the CUDA API error state so potential
312             // errors during sanity checks don't propagate.
313             if ((stat = cudaGetLastError()) != cudaSuccess)
314             {
315                 gmx_warning("An error occurred while sanity checking device #%d; %s: %s",
316                             deviceInfoList[i]->id, cudaGetErrorName(stat), cudaGetErrorString(stat));
317             }
318         }
319     }
320
321     stat = cudaPeekAtLastError();
322     GMX_RELEASE_ASSERT(stat == cudaSuccess,
323                        gmx::formatString("We promise to return with clean CUDA state, but "
324                                          "non-success state encountered: %s: %s",
325                                          cudaGetErrorName(stat), cudaGetErrorString(stat))
326                                .c_str());
327
328     return deviceInfoList;
329 }
330
331 void setActiveDevice(const DeviceInformation& deviceInfo)
332 {
333     int         deviceId = deviceInfo.id;
334     cudaError_t stat;
335
336     stat = cudaSetDevice(deviceId);
337     if (stat != cudaSuccess)
338     {
339         auto message = gmx::formatString("Failed to initialize GPU #%d", deviceId);
340         CU_RET_ERR(stat, message.c_str());
341     }
342
343     if (debug)
344     {
345         fprintf(stderr, "Initialized GPU ID #%d: %s\n", deviceId, deviceInfo.prop.name);
346     }
347 }
348
349 void releaseDevice(DeviceInformation* deviceInfo)
350 {
351     // device was used is that deviceInfo will be non-null.
352     if (deviceInfo != nullptr)
353     {
354         cudaError_t stat;
355
356         int gpuid;
357         stat = cudaGetDevice(&gpuid);
358         if (stat == cudaSuccess)
359         {
360             if (debug)
361             {
362                 fprintf(stderr, "Cleaning up context on GPU ID #%d\n", gpuid);
363             }
364
365             stat = cudaDeviceReset();
366             if (stat != cudaSuccess)
367             {
368                 gmx_warning("Failed to free GPU #%d: %s", gpuid, cudaGetErrorString(stat));
369             }
370         }
371     }
372 }
373
374 std::string getDeviceInformationString(const DeviceInformation& deviceInfo)
375 {
376     bool gpuExists = (deviceInfo.status != DeviceStatus::Nonexistent
377                       && deviceInfo.status != DeviceStatus::NonFunctional);
378
379     if (!gpuExists)
380     {
381         return gmx::formatString("#%d: %s, stat: %s", deviceInfo.id, "N/A",
382                                  c_deviceStateString[deviceInfo.status]);
383     }
384     else
385     {
386         return gmx::formatString("#%d: NVIDIA %s, compute cap.: %d.%d, ECC: %3s, stat: %s",
387                                  deviceInfo.id, deviceInfo.prop.name, deviceInfo.prop.major,
388                                  deviceInfo.prop.minor, deviceInfo.prop.ECCEnabled ? "yes" : " no",
389                                  c_deviceStateString[deviceInfo.status]);
390     }
391 }