From e30c08944cadae5b47cc31c7647e47ed136ec3b5 Mon Sep 17 00:00:00 2001 From: Artem Zhmurov Date: Wed, 24 Jun 2020 10:26:11 +0000 Subject: [PATCH] Fix minor nit-picks 1. Make the array with device status names static. 2. Expilcitely list integers in enumeration. --- src/gromacs/gpu_utils/gpu_utils.cpp | 10 +-- src/gromacs/gpu_utils/gpu_utils.cu | 82 ++++++++++--------------- src/gromacs/gpu_utils/gpu_utils.h | 3 +- src/gromacs/gpu_utils/gpu_utils_ocl.cpp | 44 ++++++------- src/gromacs/gpu_utils/gputraits.cuh | 7 ++- src/gromacs/gpu_utils/gputraits_ocl.h | 3 +- src/gromacs/hardware/CMakeLists.txt | 3 +- src/gromacs/hardware/gpu_hw_info.cpp | 52 ---------------- src/gromacs/hardware/gpu_hw_info.h | 48 +++++++++------ 9 files changed, 98 insertions(+), 154 deletions(-) delete mode 100644 src/gromacs/hardware/gpu_hw_info.cpp diff --git a/src/gromacs/gpu_utils/gpu_utils.cpp b/src/gromacs/gpu_utils/gpu_utils.cpp index 98b701ac62..fcfb1cc374 100644 --- a/src/gromacs/gpu_utils/gpu_utils.cpp +++ b/src/gromacs/gpu_utils/gpu_utils.cpp @@ -70,9 +70,9 @@ bool canPerformGpuDetection() } #if GMX_GPU == GMX_GPU_NONE -int gpu_info_get_stat(const gmx_gpu_info_t& /*unused*/, int /*unused*/) +DeviceStatus gpu_info_get_stat(const gmx_gpu_info_t& /*unused*/, int /*unused*/) { - return egpuNonexistent; + return DeviceStatus::Nonexistent; } #endif @@ -89,7 +89,7 @@ std::vector getCompatibleGpus(const gmx_gpu_info_t& gpu_info) for (int i = 0; i < gpu_info.n_dev; i++) { assert(gpu_info.deviceInfo); - if (gpu_info_get_stat(gpu_info, i) == egpuCompatible) + if (gpu_info_get_stat(gpu_info, i) == DeviceStatus::Compatible) { compatibleGpus.push_back(i); } @@ -99,8 +99,8 @@ std::vector getCompatibleGpus(const gmx_gpu_info_t& gpu_info) const char* getGpuCompatibilityDescription(const gmx_gpu_info_t& gpu_info, int index) { - return (index >= gpu_info.n_dev ? gpu_detect_res_str[egpuNonexistent] - : gpu_detect_res_str[gpu_info_get_stat(gpu_info, index)]); + return (index >= gpu_info.n_dev ? c_deviceStateString[DeviceStatus::Nonexistent] + : c_deviceStateString[gpu_info_get_stat(gpu_info, index)]); } /*! \brief Help build a descriptive message in \c error if there are * \c errorReasons why nonbondeds on a GPU are not supported. diff --git a/src/gromacs/gpu_utils/gpu_utils.cu b/src/gromacs/gpu_utils/gpu_utils.cu index b5c16c46e8..716d14588b 100644 --- a/src/gromacs/gpu_utils/gpu_utils.cu +++ b/src/gromacs/gpu_utils/gpu_utils.cu @@ -126,13 +126,19 @@ bool isHostMemoryPinned(const void* h_ptr) * Runs a series of checks to determine that the given GPU and underlying CUDA * driver/runtime functions properly. * + * \todo Currently we do not make a distinction between the type of errors + * that can appear during functionality checks. This needs to be improved, + * e.g if the dummy test kernel fails to execute with a "device busy message" + * we should appropriately report that the device is busy instead of NonFunctional. + * + * \todo Introduce errors codes and handle errors more smoothly. + * + * * \param[in] dev_id the device ID of the GPU or -1 if the device has already been initialized * \param[in] dev_prop The device properties structure * \returns 0 if the device looks OK, -1 if it sanity checks failed, and -2 if the device is busy - * - * TODO: introduce errors codes and handle errors more smoothly. */ -static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) +static DeviceStatus isDeviceFunctional(int dev_id, const cudaDeviceProp& dev_prop) { cudaError_t cu_err; int dev_count, id; @@ -141,19 +147,19 @@ static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) if (cu_err != cudaSuccess) { fprintf(stderr, "Error %d while querying device count: %s\n", cu_err, cudaGetErrorString(cu_err)); - return -1; + return DeviceStatus::NonFunctional; } /* no CUDA compatible device at all */ if (dev_count == 0) { - return -1; + return DeviceStatus::NonFunctional; } /* things might go horribly wrong if cudart is not compatible with the driver */ if (dev_count < 0 || dev_count > cuda_max_device_count) { - return -1; + return DeviceStatus::NonFunctional; } if (dev_id == -1) /* device already selected let's not destroy the context */ @@ -162,7 +168,7 @@ static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) if (cu_err != cudaSuccess) { fprintf(stderr, "Error %d while querying device id: %s\n", cu_err, cudaGetErrorString(cu_err)); - return -1; + return DeviceStatus::NonFunctional; } } else @@ -173,19 +179,19 @@ static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) fprintf(stderr, "The requested device with id %d does not seem to exist (device count=%d)\n", dev_id, dev_count); - return -1; + return DeviceStatus::NonFunctional; } } /* both major & minor is 9999 if no CUDA capable devices are present */ if (dev_prop.major == 9999 && dev_prop.minor == 9999) { - return -1; + return DeviceStatus::NonFunctional; } /* we don't care about emulation mode */ if (dev_prop.major == 0) { - return -1; + return DeviceStatus::NonFunctional; } if (id != -1) @@ -195,7 +201,7 @@ static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) { fprintf(stderr, "Error %d while switching to device #%d: %s\n", cu_err, id, cudaGetErrorString(cu_err)); - return -1; + return DeviceStatus::NonFunctional; } } @@ -205,11 +211,11 @@ static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) // if we encounter it that will happen in cudaFuncGetAttributes in the above function. if (cu_err == cudaErrorDevicesUnavailable) { - return -2; + return DeviceStatus::Unavailable; } else if (cu_err != cudaSuccess) { - return -1; + return DeviceStatus::NonFunctional; } /* try to execute a dummy kernel */ @@ -229,12 +235,12 @@ static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) fprintf(stderr, "Error occurred while running dummy kernel sanity check on device #%d:\n %s\n", id, formatExceptionMessageToString(ex).c_str()); - return -1; + return DeviceStatus::NonFunctional; } if (cudaDeviceSynchronize() != cudaSuccess) { - return -1; + return DeviceStatus::NonFunctional; } /* destroy context if we created one */ @@ -244,7 +250,7 @@ static int do_sanity_checks(int dev_id, const cudaDeviceProp& dev_prop) CU_RET_ERR(cu_err, "cudaDeviceReset failed"); } - return 0; + return DeviceStatus::Compatible; } void init_gpu(const DeviceInformation* deviceInfo) @@ -328,28 +334,13 @@ static bool is_gmx_supported_gpu(const cudaDeviceProp& dev_prop) * \param[in] deviceProp the CUDA device properties of the device checked. * \returns the status of the requested device */ -static int is_gmx_supported_gpu_id(int deviceId, const cudaDeviceProp& deviceProp) +static DeviceStatus checkDeviceStatus(int deviceId, const cudaDeviceProp& deviceProp) { if (!is_gmx_supported_gpu(deviceProp)) { - return egpuIncompatible; - } - - /* TODO: currently we do not make a distinction between the type of errors - * that can appear during sanity checks. This needs to be improved, e.g if - * the dummy test kernel fails to execute with a "device busy message" we - * should appropriately report that the device is busy instead of insane. - */ - const int checkResult = do_sanity_checks(deviceId, deviceProp); - switch (checkResult) - { - case 0: return egpuCompatible; - case -1: return egpuInsane; - case -2: return egpuUnavailable; - default: - GMX_RELEASE_ASSERT(false, "Invalid do_sanity_checks() return value"); - return egpuCompatible; + return DeviceStatus::Incompatible; } + return isDeviceFunctional(deviceId, deviceProp); } bool isGpuDetectionFunctional(std::string* errorMessage) @@ -432,22 +423,14 @@ void findGpus(gmx_gpu_info_t* gpu_info) cudaDeviceProp prop; memset(&prop, 0, sizeof(cudaDeviceProp)); stat = cudaGetDeviceProperties(&prop, i); - int checkResult; - if (stat != cudaSuccess) - { - // Will handle the error reporting below - checkResult = egpuInsane; - } - else - { - checkResult = is_gmx_supported_gpu_id(i, prop); - } + const DeviceStatus checkResult = + (stat != cudaSuccess) ? DeviceStatus::NonFunctional : checkDeviceStatus(i, prop); devs[i].id = i; devs[i].prop = prop; devs[i].stat = checkResult; - if (checkResult == egpuCompatible) + if (checkResult == DeviceStatus::Compatible) { gpu_info->n_dev_compatible++; } @@ -494,17 +477,18 @@ void get_gpu_device_info_string(char* s, const gmx_gpu_info_t& gpu_info, int ind DeviceInformation* dinfo = &gpu_info.deviceInfo[index]; - bool bGpuExists = (dinfo->stat != egpuNonexistent && dinfo->stat != egpuInsane); + bool bGpuExists = + (dinfo->stat != DeviceStatus::Nonexistent && dinfo->stat != DeviceStatus::NonFunctional); if (!bGpuExists) { - sprintf(s, "#%d: %s, stat: %s", dinfo->id, "N/A", gpu_detect_res_str[dinfo->stat]); + sprintf(s, "#%d: %s, stat: %s", dinfo->id, "N/A", c_deviceStateString[dinfo->stat]); } else { sprintf(s, "#%d: NVIDIA %s, compute cap.: %d.%d, ECC: %3s, stat: %s", dinfo->id, dinfo->prop.name, dinfo->prop.major, dinfo->prop.minor, - dinfo->prop.ECCEnabled ? "yes" : " no", gpu_detect_res_str[dinfo->stat]); + dinfo->prop.ECCEnabled ? "yes" : " no", c_deviceStateString[dinfo->stat]); } } @@ -564,7 +548,7 @@ void resetGpuProfiler(void) } } -int gpu_info_get_stat(const gmx_gpu_info_t& info, int index) +DeviceStatus gpu_info_get_stat(const gmx_gpu_info_t& info, int index) { return info.deviceInfo[index].stat; } diff --git a/src/gromacs/gpu_utils/gpu_utils.h b/src/gromacs/gpu_utils/gpu_utils.h index f935b0e513..0e27565c51 100644 --- a/src/gromacs/gpu_utils/gpu_utils.h +++ b/src/gromacs/gpu_utils/gpu_utils.h @@ -55,6 +55,7 @@ #include "gromacs/utility/basedefinitions.h" struct DeviceInformation; +enum class DeviceStatus : int; struct gmx_gpu_info_t; namespace gmx @@ -224,7 +225,7 @@ GPU_FUNC_QUALIFIER size_t sizeof_gpu_dev_info() GPU_FUNC_TERM_WITH_RETURN(0); //! Get status of device with specified index -int gpu_info_get_stat(const gmx_gpu_info_t& info, int index); +DeviceStatus gpu_info_get_stat(const gmx_gpu_info_t& info, int index); /*! \brief Check if GROMACS has been built with GPU support. * diff --git a/src/gromacs/gpu_utils/gpu_utils_ocl.cpp b/src/gromacs/gpu_utils/gpu_utils_ocl.cpp index d4e9daddd0..9f78bbda6c 100644 --- a/src/gromacs/gpu_utils/gpu_utils_ocl.cpp +++ b/src/gromacs/gpu_utils/gpu_utils_ocl.cpp @@ -140,7 +140,7 @@ static std::string makeOpenClInternalErrorString(const char* message, cl_int sta * \throws std::bad_alloc When out of memory. * \returns Whether the device passed sanity checks */ -static bool isDeviceSane(const DeviceInformation* deviceInfo, std::string* errorMessage) +static bool isDeviceFunctional(const DeviceInformation* deviceInfo, std::string* errorMessage) { cl_context_properties properties[] = { CL_CONTEXT_PLATFORM, reinterpret_cast(deviceInfo->oclPlatformId), 0 @@ -206,12 +206,12 @@ static bool isDeviceSane(const DeviceInformation* deviceInfo, std::string* error * \param[in] deviceInfo The device info pointer. * \returns The result of the compatibility checks. */ -static int isDeviceSupported(const DeviceInformation* deviceInfo) +static DeviceStatus isDeviceSupported(const DeviceInformation* deviceInfo) { if (getenv("GMX_OCL_DISABLE_COMPATIBILITY_CHECK") != nullptr) { // Assume the device is compatible because checking has been disabled. - return egpuCompatible; + return DeviceStatus::Compatible; } // OpenCL device version check, ensure >= REQUIRED_OPENCL_MIN_VERSION @@ -230,18 +230,19 @@ static int isDeviceSupported(const DeviceInformation* deviceInfo) || (deviceVersionMajor == minVersionMajor && deviceVersionMinor >= minVersionMinor))); if (!versionLargeEnough) { - return egpuIncompatible; + return DeviceStatus::Incompatible; } /* Only AMD, Intel, and NVIDIA GPUs are supported for now */ switch (deviceInfo->deviceVendor) { - case DeviceVendor::Nvidia: return egpuCompatible; + case DeviceVendor::Nvidia: return DeviceStatus::Compatible; case DeviceVendor::Amd: - return runningOnCompatibleOSForAmd() ? egpuCompatible : egpuIncompatible; + return runningOnCompatibleOSForAmd() ? DeviceStatus::Compatible : DeviceStatus::Incompatible; case DeviceVendor::Intel: - return GMX_OPENCL_NB_CLUSTER_SIZE == 4 ? egpuCompatible : egpuIncompatibleClusterSize; - default: return egpuIncompatible; + return GMX_OPENCL_NB_CLUSTER_SIZE == 4 ? DeviceStatus::Compatible + : DeviceStatus::IncompatibleClusterSize; + default: return DeviceStatus::Incompatible; } } @@ -255,26 +256,26 @@ static int isDeviceSupported(const DeviceInformation* deviceInfo) * * \param[in] deviceId The runtime-reported numeric ID of the device. * \param[in] deviceInfo The device info pointer. - * \returns An e_gpu_detect_res_t to indicate how the GPU coped with - * the sanity and compatibility check. + * \returns A DeviceStatus to indicate if the GPU device is supported and if it was able to run + * basic functionality checks. */ -static int checkGpu(size_t deviceId, const DeviceInformation* deviceInfo) +static DeviceStatus checkGpu(size_t deviceId, const DeviceInformation* deviceInfo) { - int supportStatus = isDeviceSupported(deviceInfo); - if (supportStatus != egpuCompatible) + DeviceStatus supportStatus = isDeviceSupported(deviceInfo); + if (supportStatus != DeviceStatus::Compatible) { return supportStatus; } std::string errorMessage; - if (!isDeviceSane(deviceInfo, &errorMessage)) + if (!isDeviceFunctional(deviceInfo, &errorMessage)) { gmx_warning("While sanity checking device #%zu, %s", deviceId, errorMessage.c_str()); - return egpuInsane; + return DeviceStatus::NonFunctional; } - return egpuCompatible; + return DeviceStatus::Compatible; } } // namespace gmx @@ -461,7 +462,7 @@ void findGpus(gmx_gpu_info_t* gpu_info) gpu_info->deviceInfo[device_index].stat = gmx::checkGpu(device_index, gpu_info->deviceInfo + device_index); - if (egpuCompatible == gpu_info->deviceInfo[device_index].stat) + if (DeviceStatus::Compatible == gpu_info->deviceInfo[device_index].stat) { gpu_info->n_dev_compatible++; } @@ -528,16 +529,17 @@ void get_gpu_device_info_string(char* s, const gmx_gpu_info_t& gpu_info, int ind DeviceInformation* dinfo = &gpu_info.deviceInfo[index]; - bool bGpuExists = (dinfo->stat != egpuNonexistent && dinfo->stat != egpuInsane); + bool bGpuExists = + (dinfo->stat != DeviceStatus::Nonexistent && dinfo->stat != DeviceStatus::NonFunctional); if (!bGpuExists) { - sprintf(s, "#%d: %s, stat: %s", index, "N/A", gpu_detect_res_str[dinfo->stat]); + sprintf(s, "#%d: %s, stat: %s", index, "N/A", c_deviceStateString[dinfo->stat]); } else { sprintf(s, "#%d: name: %s, vendor: %s, device version: %s, stat: %s", index, dinfo->device_name, - dinfo->vendorName, dinfo->device_version, gpu_detect_res_str[dinfo->stat]); + dinfo->vendorName, dinfo->device_version, c_deviceStateString[dinfo->stat]); } } @@ -578,7 +580,7 @@ size_t sizeof_gpu_dev_info() return sizeof(DeviceInformation); } -int gpu_info_get_stat(const gmx_gpu_info_t& info, int index) +DeviceStatus gpu_info_get_stat(const gmx_gpu_info_t& info, int index) { return info.deviceInfo[index].stat; } diff --git a/src/gromacs/gpu_utils/gputraits.cuh b/src/gromacs/gpu_utils/gputraits.cuh index 724163259d..76606611b8 100644 --- a/src/gromacs/gpu_utils/gputraits.cuh +++ b/src/gromacs/gpu_utils/gputraits.cuh @@ -44,9 +44,10 @@ * \inlibraryapi * \ingroup module_gpu_utils */ - #include +#include "gromacs/hardware/gpu_hw_info.h" + /*! \brief CUDA device information. * * The CUDA device information is queried and set at detection and contains @@ -59,8 +60,8 @@ struct DeviceInformation int id; //! CUDA device properties. cudaDeviceProp prop; - //! Result of the device check. - int stat; + //! Device status. + DeviceStatus stat; }; //! \brief Single GPU call timing event - meaningless in CUDA diff --git a/src/gromacs/gpu_utils/gputraits_ocl.h b/src/gromacs/gpu_utils/gputraits_ocl.h index c8219e5c78..ff4572e1af 100644 --- a/src/gromacs/gpu_utils/gputraits_ocl.h +++ b/src/gromacs/gpu_utils/gputraits_ocl.h @@ -46,6 +46,7 @@ */ #include "gromacs/gpu_utils/gmxopencl.h" +#include "gromacs/hardware/gpu_hw_info.h" //! OpenCL device vendors enum class DeviceVendor : int @@ -73,7 +74,7 @@ struct DeviceInformation char vendorName[256]; //!< Device vendor name. int compute_units; //!< Number of compute units. int adress_bits; //!< Number of address bits the device is capable of. - int stat; //!< Device status takes values of e_gpu_detect_res_t. + DeviceStatus stat; //!< Device status. DeviceVendor deviceVendor; //!< Device vendor. size_t maxWorkItemSizes[3]; //!< Workgroup size limits (CL_DEVICE_MAX_WORK_ITEM_SIZES). size_t maxWorkGroupSize; //!< Workgroup total size limit (CL_DEVICE_MAX_WORK_GROUP_SIZE). diff --git a/src/gromacs/hardware/CMakeLists.txt b/src/gromacs/hardware/CMakeLists.txt index 6bc268c766..efc684c31f 100644 --- a/src/gromacs/hardware/CMakeLists.txt +++ b/src/gromacs/hardware/CMakeLists.txt @@ -1,7 +1,7 @@ # # This file is part of the GROMACS molecular simulation package. # -# Copyright (c) 2015,2016,2017, by the GROMACS development team, led by +# Copyright (c) 2015,2016,2017,2020, by the GROMACS development team, led by # Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl, # and including many others, as listed in the AUTHORS file in the # top-level source directory and at http://www.gromacs.org. @@ -35,7 +35,6 @@ gmx_add_libgromacs_sources( cpuinfo.cpp detecthardware.cpp - gpu_hw_info.cpp hardwaretopology.cpp printhardware.cpp identifyavx512fmaunits.cpp diff --git a/src/gromacs/hardware/gpu_hw_info.cpp b/src/gromacs/hardware/gpu_hw_info.cpp deleted file mode 100644 index 8c79bf4320..0000000000 --- a/src/gromacs/hardware/gpu_hw_info.cpp +++ /dev/null @@ -1,52 +0,0 @@ -/* - * This file is part of the GROMACS molecular simulation package. - * - * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by - * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl, - * and including many others, as listed in the AUTHORS file in the - * top-level source directory and at http://www.gromacs.org. - * - * GROMACS is free software; you can redistribute it and/or - * modify it under the terms of the GNU Lesser General Public License - * as published by the Free Software Foundation; either version 2.1 - * of the License, or (at your option) any later version. - * - * GROMACS is distributed in the hope that it will be useful, - * but WITHOUT ANY WARRANTY; without even the implied warranty of - * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU - * Lesser General Public License for more details. - * - * You should have received a copy of the GNU Lesser General Public - * License along with GROMACS; if not, see - * http://www.gnu.org/licenses, or write to the Free Software Foundation, - * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. - * - * If you want to redistribute modifications to GROMACS, please - * consider that scientific software is very special. Version - * control is crucial - bugs must be traceable. We will be happy to - * consider code for inclusion in the official distribution, but - * derived work must not be called official GROMACS. Details are found - * in the README & COPYING files - if they are missing, get the - * official version at http://www.gromacs.org. - * - * To help us fund GROMACS development, we humbly ask that you cite - * the research papers on the package. Check out http://www.gromacs.org. - */ -#include "gmxpre.h" - -#include "gpu_hw_info.h" - -#include "config.h" - -/* Note that some of the following arrays must match the "GPU support - * enumeration" in src/config.h.cmakein, so that GMX_GPU looks up an - * array entry. */ - -// TODO If/when we unify CUDA and OpenCL support code, this should -// move to a single place in gpu_utils. -/* Names of the GPU detection/check results (see e_gpu_detect_res_t in hw_info.h). */ -const char* const gpu_detect_res_str[egpuNR] = { - "compatible", "nonexistent", - "incompatible", "incompatible (please recompile with GMX_OPENCL_NB_CLUSTER_SIZE=4)", - "insane", "unavailable" -}; diff --git a/src/gromacs/hardware/gpu_hw_info.h b/src/gromacs/hardware/gpu_hw_info.h index 8952ca4f87..aee36b4435 100644 --- a/src/gromacs/hardware/gpu_hw_info.h +++ b/src/gromacs/hardware/gpu_hw_info.h @@ -37,32 +37,40 @@ #define GMX_HARDWARE_GPU_HW_INFO_H #include "gromacs/utility/basedefinitions.h" +#include "gromacs/utility/enumerationhelpers.h" struct DeviceInformation; -/*! \brief Possible results of the GPU detection/check. - * - * The egpuInsane value means that during the sanity checks an error - * occurred that indicates malfunctioning of the device, driver, or - * incompatible driver/runtime. - * eGpuUnavailable indicates that CUDA devices are busy or unavailable - * typically due to use of cudaComputeModeExclusive, cudaComputeModeProhibited modes. - */ -typedef enum +//! Possible results of the GPU detection/check. +enum class DeviceStatus : int { - egpuCompatible = 0, - egpuNonexistent, - egpuIncompatible, - egpuIncompatibleClusterSize, - egpuInsane, - egpuUnavailable, - egpuNR -} e_gpu_detect_res_t; + //! The device is compatible + Compatible = 0, + //! Device does not exist + Nonexistent = 1, + //! Device is not compatible + Incompatible = 2, + //! OpenCL device has incompatible cluster size for non-bonded kernels. + IncompatibleClusterSize = 3, + /*! \brief An error occurred he functionality checks. + * That indicates malfunctioning of the device, driver, or incompatible driver/runtime. + */ + NonFunctional = 4, + /*! \brief CUDA devices are busy or unavailable. + * typically due to use of \p cudaComputeModeExclusive, \p cudaComputeModeProhibited modes. + */ + Unavailable = 5, + //! Enumeration size + Count = 6 +}; /*! \brief Names of the GPU detection/check results - * - * \todo Make a proper class enumeration with helper string */ -extern const char* const gpu_detect_res_str[egpuNR]; + */ +static const gmx::EnumerationArray c_deviceStateString = { + "compatible", "nonexistent", + "incompatible", "incompatible (please recompile with GMX_OPENCL_NB_CLUSTER_SIZE=4)", + "non-functional", "unavailable" +}; /*! \brief Information about GPU devices on this physical node. * -- 2.22.0