20#ifndef OPM_GPUISTL_KERNEL_ENUMS_HPP
21#define OPM_GPUISTL_KERNEL_ENUMS_HPP
23#include <opm/common/ErrorMacros.hpp>
25#include <cuda_runtime.h>
27#if CUDA_VERSION >= 12100
39 enum class MatrixStorageMPScheme {
40 DOUBLE_DIAG_DOUBLE_OFFDIAG = 0,
41 FLOAT_DIAG_FLOAT_OFFDIAG = 1,
42 DOUBLE_DIAG_FLOAT_OFFDIAG = 2
46 bool isValidMatrixStorageMPScheme(
int scheme);
49 inline MatrixStorageMPScheme makeMatrixStorageMPScheme(
int scheme) {
50 if (!detail::isValidMatrixStorageMPScheme(scheme)) {
51#if CUDA_VERSION >= 12100
52 OPM_THROW(std::invalid_argument,
53 fmt::format(fmt::runtime(
"Invalid matrix storage mixed precision scheme chosen: {}.\n"
55 "\t0: DOUBLE_DIAG_DOUBLE_OFFDIAG\n"
56 "\t1: FLOAT_DIAG_FLOAT_OFFDIAG\n"
57 "\t2: DOUBLE_DIAG_FLOAT_OFFDIAG"),
60 std::stringstream str;
61 str <<
"Invalid matrix storage mixed precision scheme chosen: " << scheme <<
".\n"
63 "\t0: DOUBLE_DIAG_DOUBLE_OFFDIAG\n"
64 "\t1: FLOAT_DIAG_FLOAT_OFFDIAG\n"
65 "\t2: DOUBLE_DIAG_FLOAT_OFFDIAG";
66 OPM_THROW(std::invalid_argument, str.str());
69 return static_cast<MatrixStorageMPScheme
>(scheme);
74 __host__ __device__
constexpr bool storeDiagonalAsFloat(MatrixStorageMPScheme scheme) {
76 case MatrixStorageMPScheme::DOUBLE_DIAG_DOUBLE_OFFDIAG:
78 case MatrixStorageMPScheme::FLOAT_DIAG_FLOAT_OFFDIAG:
80 case MatrixStorageMPScheme::DOUBLE_DIAG_FLOAT_OFFDIAG:
87 __host__ __device__
constexpr bool storeOffDiagonalAsFloat(MatrixStorageMPScheme scheme) {
89 case MatrixStorageMPScheme::DOUBLE_DIAG_DOUBLE_OFFDIAG:
91 case MatrixStorageMPScheme::FLOAT_DIAG_FLOAT_OFFDIAG:
93 case MatrixStorageMPScheme::DOUBLE_DIAG_FLOAT_OFFDIAG:
101 __host__ __device__
constexpr bool usingMixedPrecision(MatrixStorageMPScheme scheme) {
103 case MatrixStorageMPScheme::DOUBLE_DIAG_DOUBLE_OFFDIAG:
105 case MatrixStorageMPScheme::FLOAT_DIAG_FLOAT_OFFDIAG:
107 case MatrixStorageMPScheme::DOUBLE_DIAG_FLOAT_OFFDIAG:
114 inline bool isValidMatrixStorageMPScheme(
int scheme) {
115 switch (
static_cast<MatrixStorageMPScheme
>(scheme)) {
116 case MatrixStorageMPScheme::DOUBLE_DIAG_DOUBLE_OFFDIAG:
117 case MatrixStorageMPScheme::FLOAT_DIAG_FLOAT_OFFDIAG:
118 case MatrixStorageMPScheme::DOUBLE_DIAG_FLOAT_OFFDIAG:
Contains wrappers to make the CuBLAS library behave as a modern C++ library with function overlading.
Definition autotuner.hpp:30
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition AmgxInterface.hpp:38