19#ifndef OPM_GPUSPARSEMATRIXWRAPPER_HPP
20#define OPM_GPUSPARSEMATRIXWRAPPER_HPP
26#include <opm/common/ErrorMacros.hpp>
27#include <opm/simulators/linalg/gpuistl/GpuVector.hpp>
28#include <opm/simulators/linalg/gpuistl/GpuSparseMatrix.hpp>
29#include <opm/simulators/linalg/gpuistl/GpuSparseMatrixGeneric.hpp>
30#include <opm/simulators/linalg/gpuistl/detail/CuMatrixDescription.hpp>
31#include <opm/simulators/linalg/gpuistl/detail/CuSparseHandle.hpp>
32#include <opm/simulators/linalg/gpuistl/detail/safe_conversion.hpp>
71#if USE_HIP || (!USE_HIP && CUDA_VERSION < 13000)
78 matrix_type* operator->() {
80 throw std::runtime_error(
"GpuSparseMatrixWrapper: underlying matrix is nullptr.");
82 return m_matrix.get();
84 const matrix_type* operator->()
const {
86 throw std::runtime_error(
"GpuSparseMatrixWrapper: underlying matrix is nullptr.");
88 return m_matrix.get();
115 const int* rowIndices,
116 const int* columnIndices,
117 size_t numberOfNonzeroBlocks,
121 m_matrix = std::make_unique<matrix_type>(nonZeroElements,
124 numberOfNonzeroBlocks,
143 m_matrix = std::make_unique<matrix_type>(rowIndices, columnIndices,
blockSize);
148 if (!other.m_matrix) {
149 throw std::runtime_error(
"Internal error, other.m_matrix is a nullptr.");
151 m_matrix = std::make_unique<matrix_type>(*other.m_matrix);
158 ~GpuSparseMatrixWrapper() =
default;
162 const matrix_type& get()
const {
164 throw std::runtime_error(
"GpuSparseMatrixWrapper: underlying matrix is nullptr.");
178 template <
class MatrixType>
182 gpuSparseMatrixWrapper.m_matrix = std::make_unique<matrix_type>(
184 return gpuSparseMatrixWrapper;
188 template <
class M = matrix_type,
189 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
190 void setUpperTriangular()
192 m_matrix->setUpperTriangular();
198 template <
class M = matrix_type,
199 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
202 m_matrix->setLowerTriangular();
208 template <
class M = matrix_type,
209 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
212 m_matrix->setUnitDiagonal();
218 template <
class M = matrix_type,
219 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
222 m_matrix->setNonUnitDiagonal();
230 return m_matrix->N();
239 return m_matrix->nonzeroes();
249 return m_matrix->getNonZeroValues();
259 return m_matrix->getNonZeroValues();
269 return m_matrix->getRowIndices();
279 return m_matrix->getRowIndices();
289 return m_matrix->getColumnIndices();
299 return m_matrix->getColumnIndices();
310 return m_matrix->dim();
318 return m_matrix->blockSize();
328 return m_matrix->getDescription();
360 m_matrix->usmv(alpha, x, y);
373 template <
class MatrixType>
376 m_matrix->updateNonzeroValues(matrix, copyNonZeroElementsDirectly);
386 m_matrix->updateNonzeroValues(matrix.get());
410 template<
class FunctionType>
413 return m_matrix->dispatchOnBlocksize(function);
418 std::unique_ptr<matrix_type> m_matrix =
nullptr;
The GpuSparseMatrixGeneric class uses cuSPARSE Generic API for sparse matrix operations.
Definition GpuSparseMatrixGeneric.hpp:51
The GpuSparseMatrixWrapper Checks CUDA/HIP version and dispatches a version either using the old or t...
Definition GpuSparseMatrixWrapper.hpp:61
GpuVector< int > & getRowIndices()
getRowIndices returns the row indices used to represent the BSR structure.
Definition GpuSparseMatrixWrapper.hpp:267
virtual void usmv(T alpha, const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=alpha * Ax + y
Definition GpuSparseMatrixWrapper.hpp:358
auto dispatchOnBlocksize(FunctionType function) const
Dispatches a function based on the block size of the matrix.
Definition GpuSparseMatrixWrapper.hpp:411
void updateNonzeroValues(const GpuSparseMatrixWrapper< T > &matrix)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
Definition GpuSparseMatrixWrapper.hpp:384
detail::GpuSparseMatrixDescription & getDescription()
getDescription the cusparse matrix description.
Definition GpuSparseMatrixWrapper.hpp:326
static GpuSparseMatrixWrapper< T > fromMatrix(const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
fromMatrix creates a new matrix with the same block size and values as the given matrix
Definition GpuSparseMatrixWrapper.hpp:179
GpuVector< int > & getColumnIndices()
getColumnIndices returns the column indices used to represent the BSR structure.
Definition GpuSparseMatrixWrapper.hpp:287
virtual void mv(const GpuVector< T > &x, GpuVector< T > &y) const
mv performs matrix vector multiply y = Ax
Definition GpuSparseMatrixWrapper.hpp:336
void setUnitDiagonal()
setUnitDiagonal sets the CuSparse flag that this has unit diagional.
Definition GpuSparseMatrixWrapper.hpp:210
GpuVector< T > & getNonZeroValues()
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition GpuSparseMatrixWrapper.hpp:247
void setLowerTriangular()
setLowerTriangular sets the CuSparse flag that this is an lower diagonal (with non-unit diagonal) mat...
Definition GpuSparseMatrixWrapper.hpp:200
size_t blockSize() const
Definition GpuSparseMatrixWrapper.hpp:316
const GpuVector< int > & getRowIndices() const
getRowIndices returns the row indices used to represent the BSR structure.
Definition GpuSparseMatrixWrapper.hpp:277
size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition GpuSparseMatrixWrapper.hpp:237
const GpuVector< T > & getNonZeroValues() const
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition GpuSparseMatrixWrapper.hpp:257
GpuSparseMatrixWrapper(const T *nonZeroElements, const int *rowIndices, const int *columnIndices, size_t numberOfNonzeroBlocks, size_t blockSize, size_t numberOfRows)
Create the sparse matrix specified by the raw data.
Definition GpuSparseMatrixWrapper.hpp:114
const GpuVector< int > & getColumnIndices() const
getColumnIndices returns the column indices used to represent the BSR structure.
Definition GpuSparseMatrixWrapper.hpp:297
virtual void umv(const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=Ax+y
Definition GpuSparseMatrixWrapper.hpp:346
void updateNonzeroValues(const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
Definition GpuSparseMatrixWrapper.hpp:374
void setNonUnitDiagonal()
setNonUnitDiagonal sets the CuSparse flag that this has non-unit diagional.
Definition GpuSparseMatrixWrapper.hpp:220
static constexpr int max_block_size
Definition GpuSparseMatrixWrapper.hpp:99
size_t N() const
N returns the number of rows (which is equal to the number of columns).
Definition GpuSparseMatrixWrapper.hpp:228
GpuSparseMatrixWrapper(const GpuVector< int > &rowIndices, const GpuVector< int > &columnIndices, size_t blockSize)
Create a sparse matrix by copying the sparsity structure of another matrix, not filling in the values...
Definition GpuSparseMatrixWrapper.hpp:139
size_t dim() const
dim returns the dimension of the vector space on which this matrix acts
Definition GpuSparseMatrixWrapper.hpp:308
The GpuSparseMatrix class simple wrapper class for a CuSparse matrix.
Definition GpuSparseMatrix.hpp:60
static GpuSparseMatrix< T > fromMatrix(const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
fromMatrix creates a new matrix with the same block size and values as the given matrix
Definition GpuSparseMatrix.cpp:115
Definition gpu_type_detection.hpp:30
CuSparseResource< cusparseMatDescr_t > GpuSparseMatrixDescription
GpuSparseMatrixDescription holder.
Definition CuMatrixDescription.hpp:30