3cf2eec706b38dd375f384686672cdb09252bc5e
[alexxy/gromacs.git] / src / gromacs / hardware / device_management_ocl.cpp
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 Define functions for detection and initialization for OpenCL devices.
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  */
45
46 #include "gmxpre.h"
47
48 #include "config.h"
49
50 #include <assert.h>
51 #include <stdio.h>
52 #include <stdlib.h>
53 #include <string.h>
54
55 #include <cstdio>
56 #ifdef __APPLE__
57 #    include <sys/sysctl.h>
58 #endif
59
60 #include <memory.h>
61
62 #include "gromacs/gpu_utils/ocl_compiler.h"
63 #include "gromacs/gpu_utils/oclraii.h"
64 #include "gromacs/gpu_utils/oclutils.h"
65 #include "gromacs/hardware/device_information.h"
66 #include "gromacs/hardware/device_management.h"
67 #include "gromacs/hardware/hw_info.h"
68 #include "gromacs/utility/cstringutil.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/smalloc.h"
72 #include "gromacs/utility/stringutil.h"
73
74 /*! \brief Return true if executing on compatible OS for AMD OpenCL.
75  *
76  * This is assumed to be true for OS X version of at least 10.10.4 and
77  * all other OS flavors.
78  *
79  * Uses the BSD sysctl() interfaces to extract the kernel version.
80  *
81  * \return true if version is 14.4 or later (= OS X version 10.10.4),
82  *         or OS is not Darwin.
83  */
84 static bool runningOnCompatibleOSForAmd()
85 {
86 #ifdef __APPLE__
87     int    mib[2];
88     char   kernelVersion[256];
89     size_t len = sizeof(kernelVersion);
90
91     mib[0] = CTL_KERN;
92     mib[1] = KERN_OSRELEASE;
93
94     sysctl(mib, sizeof(mib) / sizeof(mib[0]), kernelVersion, &len, NULL, 0);
95
96     int major = strtod(kernelVersion, NULL);
97     int minor = strtod(strchr(kernelVersion, '.') + 1, NULL);
98
99     // Kernel 14.4 corresponds to OS X 10.10.4
100     return (major > 14 || (major == 14 && minor >= 4));
101 #else
102     return true;
103 #endif
104 }
105
106 namespace gmx
107 {
108
109 /*! \brief Make an error string following an OpenCL API call.
110  *
111  *  It is meant to be called with \p status != CL_SUCCESS, but it will
112  *  work correctly even if it is called with no OpenCL failure.
113  *
114  * \param[in]  message  Supplies context, e.g. the name of the API call that returned the error.
115  * \param[in]  status   OpenCL API status code
116  * \returns             A string describing the OpenCL error.
117  */
118 static std::string makeOpenClInternalErrorString(const char* message, cl_int status)
119 {
120     if (message != nullptr)
121     {
122         return formatString("%s did %ssucceed %d: %s", message, ((status != CL_SUCCESS) ? "not " : ""),
123                             status, ocl_get_error_string(status).c_str());
124     }
125     else
126     {
127         return formatString("%sOpenCL error encountered %d: %s", ((status != CL_SUCCESS) ? "" : "No "),
128                             status, ocl_get_error_string(status).c_str());
129     }
130 }
131
132 /*!
133  * \brief Checks that device \c deviceInfo is sane (ie can run a kernel).
134  *
135  * Compiles and runs a dummy kernel to determine whether the given
136  * OpenCL device functions properly.
137  *
138  *
139  * \param[in]  deviceInfo      The device info pointer.
140  * \param[out] errorMessage    An error message related to a failing OpenCL API call.
141  * \throws     std::bad_alloc  When out of memory.
142  * \returns                    Whether the device passed sanity checks
143  */
144 static bool isDeviceFunctional(const DeviceInformation* deviceInfo, std::string* errorMessage)
145 {
146     cl_context_properties properties[] = {
147         CL_CONTEXT_PLATFORM, reinterpret_cast<cl_context_properties>(deviceInfo->oclPlatformId), 0
148     };
149     // uncrustify spacing
150
151     cl_int    status;
152     auto      deviceId = deviceInfo->oclDeviceId;
153     ClContext context(clCreateContext(properties, 1, &deviceId, nullptr, nullptr, &status));
154     if (status != CL_SUCCESS)
155     {
156         errorMessage->assign(makeOpenClInternalErrorString("clCreateContext", status));
157         return false;
158     }
159     ClCommandQueue commandQueue(clCreateCommandQueue(context, deviceId, 0, &status));
160     if (status != CL_SUCCESS)
161     {
162         errorMessage->assign(makeOpenClInternalErrorString("clCreateCommandQueue", status));
163         return false;
164     }
165
166     // Some compilers such as Apple's require kernel functions to have at least one argument
167     const char* lines[] = { "__kernel void dummyKernel(__global void* input){}" };
168     ClProgram   program(clCreateProgramWithSource(context, 1, lines, nullptr, &status));
169     if (status != CL_SUCCESS)
170     {
171         errorMessage->assign(makeOpenClInternalErrorString("clCreateProgramWithSource", status));
172         return false;
173     }
174
175     if ((status = clBuildProgram(program, 0, nullptr, nullptr, nullptr, nullptr)) != CL_SUCCESS)
176     {
177         errorMessage->assign(makeOpenClInternalErrorString("clBuildProgram", status));
178         return false;
179     }
180
181     ClKernel kernel(clCreateKernel(program, "dummyKernel", &status));
182     if (status != CL_SUCCESS)
183     {
184         errorMessage->assign(makeOpenClInternalErrorString("clCreateKernel", status));
185         return false;
186     }
187
188     clSetKernelArg(kernel, 0, sizeof(void*), nullptr);
189
190     const size_t localWorkSize = 1, globalWorkSize = 1;
191     if ((status = clEnqueueNDRangeKernel(commandQueue, kernel, 1, nullptr, &globalWorkSize,
192                                          &localWorkSize, 0, nullptr, nullptr))
193         != CL_SUCCESS)
194     {
195         errorMessage->assign(makeOpenClInternalErrorString("clEnqueueNDRangeKernel", status));
196         return false;
197     }
198     return true;
199 }
200
201 /*!
202  * \brief Checks that device \c deviceInfo is compatible with GROMACS.
203  *
204  *  Vendor and OpenCL version support checks are executed an the result
205  *  of these returned.
206  *
207  * \param[in]  deviceInfo  The device info pointer.
208  * \returns                The result of the compatibility checks.
209  */
210 static DeviceStatus isDeviceSupported(const DeviceInformation* deviceInfo)
211 {
212     if (getenv("GMX_OCL_DISABLE_COMPATIBILITY_CHECK") != nullptr)
213     {
214         // Assume the device is compatible because checking has been disabled.
215         return DeviceStatus::Compatible;
216     }
217
218     // OpenCL device version check, ensure >= REQUIRED_OPENCL_MIN_VERSION
219     constexpr unsigned int minVersionMajor = REQUIRED_OPENCL_MIN_VERSION_MAJOR;
220     constexpr unsigned int minVersionMinor = REQUIRED_OPENCL_MIN_VERSION_MINOR;
221
222     // Based on the OpenCL spec we're checking the version supported by
223     // the device which has the following format:
224     //      OpenCL<space><major_version.minor_version><space><vendor-specific information>
225     unsigned int deviceVersionMinor, deviceVersionMajor;
226     const int    valuesScanned = std::sscanf(deviceInfo->device_version, "OpenCL %u.%u",
227                                           &deviceVersionMajor, &deviceVersionMinor);
228     const bool   versionLargeEnough =
229             ((valuesScanned == 2)
230              && ((deviceVersionMajor > minVersionMajor)
231                  || (deviceVersionMajor == minVersionMajor && deviceVersionMinor >= minVersionMinor)));
232     if (!versionLargeEnough)
233     {
234         return DeviceStatus::Incompatible;
235     }
236
237     /* Only AMD, Intel, and NVIDIA GPUs are supported for now */
238     switch (deviceInfo->deviceVendor)
239     {
240         case DeviceVendor::Nvidia: return DeviceStatus::Compatible;
241         case DeviceVendor::Amd:
242             return runningOnCompatibleOSForAmd() ? DeviceStatus::Compatible : DeviceStatus::Incompatible;
243         case DeviceVendor::Intel:
244             return GMX_OPENCL_NB_CLUSTER_SIZE == 4 ? DeviceStatus::Compatible
245                                                    : DeviceStatus::IncompatibleClusterSize;
246         default: return DeviceStatus::Incompatible;
247     }
248 }
249
250
251 /*! \brief Check whether the \c ocl_gpu_device is suitable for use by mdrun
252  *
253  * Runs sanity checks: checking that the runtime can compile a dummy kernel
254  * and this can be executed;
255  * Runs compatibility checks verifying the device OpenCL version requirement
256  * and vendor/OS support.
257  *
258  * \param[in]  deviceId      The runtime-reported numeric ID of the device.
259  * \param[in]  deviceInfo    The device info pointer.
260  * \returns  A DeviceStatus to indicate if the GPU device is supported and if it was able to run
261  *           basic functionality checks.
262  */
263 static DeviceStatus checkGpu(size_t deviceId, const DeviceInformation* deviceInfo)
264 {
265
266     DeviceStatus supportStatus = isDeviceSupported(deviceInfo);
267     if (supportStatus != DeviceStatus::Compatible)
268     {
269         return supportStatus;
270     }
271
272     std::string errorMessage;
273     if (!isDeviceFunctional(deviceInfo, &errorMessage))
274     {
275         gmx_warning("While sanity checking device #%zu, %s", deviceId, errorMessage.c_str());
276         return DeviceStatus::NonFunctional;
277     }
278
279     return DeviceStatus::Compatible;
280 }
281
282 } // namespace gmx
283
284 /*! \brief Returns an DeviceVendor value corresponding to the input OpenCL vendor name.
285  *
286  *  \param[in] vendorName  String with OpenCL vendor name.
287  *  \returns               DeviceVendor value for the input vendor name
288  */
289 static DeviceVendor getDeviceVendor(const char* vendorName)
290 {
291     if (vendorName)
292     {
293         if (strstr(vendorName, "NVIDIA"))
294         {
295             return DeviceVendor::Nvidia;
296         }
297         else if (strstr(vendorName, "AMD") || strstr(vendorName, "Advanced Micro Devices"))
298         {
299             return DeviceVendor::Amd;
300         }
301         else if (strstr(vendorName, "Intel"))
302         {
303             return DeviceVendor::Intel;
304         }
305     }
306     return DeviceVendor::Unknown;
307 }
308
309 bool isGpuDetectionFunctional(std::string* errorMessage)
310 {
311     cl_uint numPlatforms;
312     cl_int  status = clGetPlatformIDs(0, nullptr, &numPlatforms);
313     GMX_ASSERT(status != CL_INVALID_VALUE, "Incorrect call of clGetPlatformIDs detected");
314 #ifdef cl_khr_icd
315     if (status == CL_PLATFORM_NOT_FOUND_KHR)
316     {
317         // No valid ICDs found
318         if (errorMessage != nullptr)
319         {
320             errorMessage->assign("No valid OpenCL driver found");
321         }
322         return false;
323     }
324 #endif
325     GMX_RELEASE_ASSERT(
326             status == CL_SUCCESS,
327             gmx::formatString("An unexpected value was returned from clGetPlatformIDs %d: %s",
328                               status, ocl_get_error_string(status).c_str())
329                     .c_str());
330     bool foundPlatform = (numPlatforms > 0);
331     if (!foundPlatform && errorMessage != nullptr)
332     {
333         errorMessage->assign("No OpenCL platforms found even though the driver was valid");
334     }
335     return foundPlatform;
336 }
337
338 void findGpus(gmx_gpu_info_t* gpu_info)
339 {
340     cl_uint         ocl_platform_count;
341     cl_platform_id* ocl_platform_ids;
342     cl_device_type  req_dev_type = CL_DEVICE_TYPE_GPU;
343
344     ocl_platform_ids = nullptr;
345
346     if (getenv("GMX_OCL_FORCE_CPU") != nullptr)
347     {
348         req_dev_type = CL_DEVICE_TYPE_CPU;
349     }
350
351     while (true)
352     {
353         cl_int status = clGetPlatformIDs(0, nullptr, &ocl_platform_count);
354         if (CL_SUCCESS != status)
355         {
356             GMX_THROW(gmx::InternalError(
357                     gmx::formatString("An unexpected value %d was returned from clGetPlatformIDs: ", status)
358                     + ocl_get_error_string(status)));
359         }
360
361         if (1 > ocl_platform_count)
362         {
363             // TODO this should have a descriptive error message that we only support one OpenCL platform
364             break;
365         }
366
367         snew(ocl_platform_ids, ocl_platform_count);
368
369         status = clGetPlatformIDs(ocl_platform_count, ocl_platform_ids, nullptr);
370         if (CL_SUCCESS != status)
371         {
372             GMX_THROW(gmx::InternalError(
373                     gmx::formatString("An unexpected value %d was returned from clGetPlatformIDs: ", status)
374                     + ocl_get_error_string(status)));
375         }
376
377         for (unsigned int i = 0; i < ocl_platform_count; i++)
378         {
379             cl_uint ocl_device_count;
380
381             /* If requesting req_dev_type devices fails, just go to the next platform */
382             if (CL_SUCCESS != clGetDeviceIDs(ocl_platform_ids[i], req_dev_type, 0, nullptr, &ocl_device_count))
383             {
384                 continue;
385             }
386
387             if (1 <= ocl_device_count)
388             {
389                 gpu_info->n_dev += ocl_device_count;
390             }
391         }
392
393         if (1 > gpu_info->n_dev)
394         {
395             break;
396         }
397
398         snew(gpu_info->deviceInfo, gpu_info->n_dev);
399
400         {
401             int           device_index;
402             cl_device_id* ocl_device_ids;
403
404             snew(ocl_device_ids, gpu_info->n_dev);
405             device_index = 0;
406
407             for (unsigned int i = 0; i < ocl_platform_count; i++)
408             {
409                 cl_uint ocl_device_count;
410
411                 /* If requesting req_dev_type devices fails, just go to the next platform */
412                 if (CL_SUCCESS
413                     != clGetDeviceIDs(ocl_platform_ids[i], req_dev_type, gpu_info->n_dev,
414                                       ocl_device_ids, &ocl_device_count))
415                 {
416                     continue;
417                 }
418
419                 if (1 > ocl_device_count)
420                 {
421                     break;
422                 }
423
424                 for (unsigned int j = 0; j < ocl_device_count; j++)
425                 {
426                     gpu_info->deviceInfo[device_index].oclPlatformId = ocl_platform_ids[i];
427                     gpu_info->deviceInfo[device_index].oclDeviceId   = ocl_device_ids[j];
428
429                     gpu_info->deviceInfo[device_index].device_name[0] = 0;
430                     clGetDeviceInfo(ocl_device_ids[j], CL_DEVICE_NAME,
431                                     sizeof(gpu_info->deviceInfo[device_index].device_name),
432                                     gpu_info->deviceInfo[device_index].device_name, nullptr);
433
434                     gpu_info->deviceInfo[device_index].device_version[0] = 0;
435                     clGetDeviceInfo(ocl_device_ids[j], CL_DEVICE_VERSION,
436                                     sizeof(gpu_info->deviceInfo[device_index].device_version),
437                                     gpu_info->deviceInfo[device_index].device_version, nullptr);
438
439                     gpu_info->deviceInfo[device_index].vendorName[0] = 0;
440                     clGetDeviceInfo(ocl_device_ids[j], CL_DEVICE_VENDOR,
441                                     sizeof(gpu_info->deviceInfo[device_index].vendorName),
442                                     gpu_info->deviceInfo[device_index].vendorName, nullptr);
443
444                     gpu_info->deviceInfo[device_index].compute_units = 0;
445                     clGetDeviceInfo(ocl_device_ids[j], CL_DEVICE_MAX_COMPUTE_UNITS,
446                                     sizeof(gpu_info->deviceInfo[device_index].compute_units),
447                                     &(gpu_info->deviceInfo[device_index].compute_units), nullptr);
448
449                     gpu_info->deviceInfo[device_index].adress_bits = 0;
450                     clGetDeviceInfo(ocl_device_ids[j], CL_DEVICE_ADDRESS_BITS,
451                                     sizeof(gpu_info->deviceInfo[device_index].adress_bits),
452                                     &(gpu_info->deviceInfo[device_index].adress_bits), nullptr);
453
454                     gpu_info->deviceInfo[device_index].deviceVendor =
455                             getDeviceVendor(gpu_info->deviceInfo[device_index].vendorName);
456
457                     clGetDeviceInfo(ocl_device_ids[j], CL_DEVICE_MAX_WORK_ITEM_SIZES, 3 * sizeof(size_t),
458                                     &gpu_info->deviceInfo[device_index].maxWorkItemSizes, nullptr);
459
460                     clGetDeviceInfo(ocl_device_ids[j], CL_DEVICE_MAX_WORK_GROUP_SIZE, sizeof(size_t),
461                                     &gpu_info->deviceInfo[device_index].maxWorkGroupSize, nullptr);
462
463                     gpu_info->deviceInfo[device_index].stat =
464                             gmx::checkGpu(device_index, gpu_info->deviceInfo + device_index);
465
466                     if (DeviceStatus::Compatible == gpu_info->deviceInfo[device_index].stat)
467                     {
468                         gpu_info->n_dev_compatible++;
469                     }
470
471                     device_index++;
472                 }
473             }
474
475             gpu_info->n_dev = device_index;
476
477             /* Dummy sort of devices -  AMD first, then NVIDIA, then Intel */
478             // TODO: Sort devices based on performance.
479             if (0 < gpu_info->n_dev)
480             {
481                 int last = -1;
482                 for (int i = 0; i < gpu_info->n_dev; i++)
483                 {
484                     if (gpu_info->deviceInfo[i].deviceVendor == DeviceVendor::Amd)
485                     {
486                         last++;
487
488                         if (last < i)
489                         {
490                             std::swap(gpu_info->deviceInfo[i], gpu_info->deviceInfo[last]);
491                         }
492                     }
493                 }
494
495                 /* if more than 1 device left to be sorted */
496                 if ((gpu_info->n_dev - 1 - last) > 1)
497                 {
498                     for (int i = 0; i < gpu_info->n_dev; i++)
499                     {
500                         if (gpu_info->deviceInfo[i].deviceVendor == DeviceVendor::Nvidia)
501                         {
502                             last++;
503
504                             if (last < i)
505                             {
506                                 std::swap(gpu_info->deviceInfo[i], gpu_info->deviceInfo[last]);
507                             }
508                         }
509                     }
510                 }
511             }
512
513             sfree(ocl_device_ids);
514         }
515
516         break;
517     }
518
519     sfree(ocl_platform_ids);
520 }
521
522 void init_gpu(const DeviceInformation* deviceInfo)
523 {
524     assert(deviceInfo);
525
526     // If the device is NVIDIA, for safety reasons we disable the JIT
527     // caching as this is known to be broken at least until driver 364.19;
528     // the cache does not always get regenerated when the source code changes,
529     // e.g. if the path to the kernel sources remains the same
530
531     if (deviceInfo->deviceVendor == DeviceVendor::Nvidia)
532     {
533         // Ignore return values, failing to set the variable does not mean
534         // that something will go wrong later.
535 #ifdef _MSC_VER
536         _putenv("CUDA_CACHE_DISABLE=1");
537 #else
538         // Don't override, maybe a dev is testing.
539         setenv("CUDA_CACHE_DISABLE", "1", 0);
540 #endif
541     }
542 }
543
544 void free_gpu(const DeviceInformation* /* deviceInfo */) {}
545
546 DeviceInformation* getDeviceInfo(const gmx_gpu_info_t& gpu_info, int deviceId)
547 {
548     if (deviceId < 0 || deviceId >= gpu_info.n_dev)
549     {
550         gmx_incons("Invalid GPU deviceId requested");
551     }
552     return &gpu_info.deviceInfo[deviceId];
553 }
554
555 void get_gpu_device_info_string(char* s, const gmx_gpu_info_t& gpu_info, int index)
556 {
557     assert(s);
558
559     if (index < 0 && index >= gpu_info.n_dev)
560     {
561         return;
562     }
563
564     DeviceInformation* dinfo = &gpu_info.deviceInfo[index];
565
566     bool bGpuExists =
567             (dinfo->stat != DeviceStatus::Nonexistent && dinfo->stat != DeviceStatus::NonFunctional);
568
569     if (!bGpuExists)
570     {
571         sprintf(s, "#%d: %s, stat: %s", index, "N/A", c_deviceStateString[dinfo->stat]);
572     }
573     else
574     {
575         sprintf(s, "#%d: name: %s, vendor: %s, device version: %s, stat: %s", index, dinfo->device_name,
576                 dinfo->vendorName, dinfo->device_version, c_deviceStateString[dinfo->stat]);
577     }
578 }
579
580 size_t sizeof_gpu_dev_info()
581 {
582     return sizeof(DeviceInformation);
583 }
584
585 DeviceStatus gpu_info_get_stat(const gmx_gpu_info_t& info, int index)
586 {
587     return info.deviceInfo[index].stat;
588 }