GpuSparseMatrixWrapper.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 Equinor ASA
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef OPM_GPUSPARSEMATRIXWRAPPER_HPP
20#define OPM_GPUSPARSEMATRIXWRAPPER_HPP
21
22#include <opm/common/ErrorMacros.hpp>
29
30#include <cstddef>
31#include <cuda.h>
32#include <cusparse.h>
33#include <memory>
34#include <stdexcept>
35#include <type_traits>
36
37namespace Opm::gpuistl
38{
39
60template <typename T, bool ForceLegacy = false>
62{
63public:
64 using field_type = T;
65
66 /*
67 Here is the primary function of this class.
68 Since the generic API for CUDA/HIP is primaryly supported on CUDA 13 (and not yet HIP) for blocked
69 matrices, places wanting to use blocked matrices can invoke this class which handles which API to use.
70 Basically we just check if HIP is present, or if we are using cuda and a version prior to 13.
71 If ForceLegacy is true, always use the legacy API (GpuSparseMatrix) regardless of CUDA version.
72 */
73#if USE_HIP || (!USE_HIP && CUDA_VERSION < 13000)
75#else
76 using matrix_type = std::conditional_t<ForceLegacy,
79#endif
80
81 // Arrow operator overloads for direct access to the underlying matrix
83 if (!m_matrix) {
84 throw std::runtime_error("GpuSparseMatrixWrapper: underlying matrix is nullptr.");
85 }
86 return m_matrix.get();
87 }
88 const matrix_type* operator->() const {
89 if (!m_matrix) {
90 throw std::runtime_error("GpuSparseMatrixWrapper: underlying matrix is nullptr.");
91 }
92 return m_matrix.get();
93 }
94
103 static constexpr int max_block_size = 6;
104
118 GpuSparseMatrixWrapper(const T* nonZeroElements,
119 const int* rowIndices,
120 const int* columnIndices,
121 std::size_t numberOfNonzeroBlocks,
122 std::size_t blockSize,
123 std::size_t numberOfRows)
124 {
125 m_matrix = std::make_unique<matrix_type>(nonZeroElements,
126 rowIndices,
127 columnIndices,
128 numberOfNonzeroBlocks,
129 blockSize,
130 numberOfRows);
131 }
132
144 const GpuVector<int>& columnIndices,
145 std::size_t blockSize)
146 {
147 m_matrix = std::make_unique<matrix_type>(rowIndices, columnIndices, blockSize);
148 }
149
151 {
152 if (!other.m_matrix) {
153 throw std::runtime_error("Internal error, other.m_matrix is a nullptr.");
154 }
155 m_matrix = std::make_unique<matrix_type>(*other.m_matrix);
156 }
157
158 // We want to have this as non-mutable as possible, that is we do not want
159 // to deal with changing matrix sizes and sparsity patterns.
161
162 virtual ~GpuSparseMatrixWrapper() = default;
163
165
166 const matrix_type& get() const {
167 if (!m_matrix) {
168 throw std::runtime_error("GpuSparseMatrixWrapper: underlying matrix is nullptr.");
169 }
170 return *m_matrix;
171 }
172
182 template <class MatrixType>
183 static GpuSparseMatrixWrapper<T, ForceLegacy> fromMatrix(const MatrixType& matrix, bool copyNonZeroElementsDirectly = false)
184 {
185 GpuSparseMatrixWrapper<T, ForceLegacy> gpuSparseMatrixWrapper;
186 gpuSparseMatrixWrapper.m_matrix = std::make_unique<matrix_type>(
187 matrix_type::fromMatrix(matrix, copyNonZeroElementsDirectly));
188 return gpuSparseMatrixWrapper;
189 }
190
191 // Only participates in overload resolution when matrix_type == GpuSparseMatrix<T>
192 template <class M = matrix_type,
193 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
195 {
196 m_matrix->setUpperTriangular();
197 }
198
202 template <class M = matrix_type,
203 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
205 {
206 m_matrix->setLowerTriangular();
207 }
208
212 template <class M = matrix_type,
213 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
215 {
216 m_matrix->setUnitDiagonal();
217 }
218
222 template <class M = matrix_type,
223 typename = std::enable_if_t<std::is_same_v<M, GpuSparseMatrix<T>>>>
225 {
226 m_matrix->setNonUnitDiagonal();
227 }
228
232 std::size_t N() const
233 {
234 return m_matrix->N();
235 }
236
241 std::size_t nonzeroes() const
242 {
243 return m_matrix->nonzeroes();
244 }
245
252 {
253 return m_matrix->getNonZeroValues();
254 }
255
262 {
263 return m_matrix->getNonZeroValues();
264 }
265
272 {
273 return m_matrix->getRowIndices();
274 }
275
282 {
283 return m_matrix->getRowIndices();
284 }
285
292 {
293 return m_matrix->getColumnIndices();
294 }
295
302 {
303 return m_matrix->getColumnIndices();
304 }
305
312 std::size_t dim() const
313 {
314 return m_matrix->dim();
315 }
316
320 std::size_t blockSize() const
321 {
322 return m_matrix->blockSize();
323 }
324
331 {
332 return m_matrix->getDescription();
333 }
334
340 virtual void mv(const GpuVector<T>& x, GpuVector<T>& y) const
341 {
342 m_matrix->mv(x, y);
343 }
344
350 virtual void umv(const GpuVector<T>& x, GpuVector<T>& y) const
351 {
352 m_matrix->umv(x, y);
353 }
354
355
362 virtual void usmv(T alpha, const GpuVector<T>& x, GpuVector<T>& y) const
363 {
364 m_matrix->usmv(alpha, x, y);
365 }
366
377 template <class MatrixType>
378 void updateNonzeroValues(const MatrixType& matrix, bool copyNonZeroElementsDirectly = false)
379 {
380 m_matrix->updateNonzeroValues(matrix, copyNonZeroElementsDirectly);
381 }
382
388 template <bool OtherForceLegacy>
390 {
391 m_matrix->updateNonzeroValues(matrix.get());
392 }
393
394
415 template<class FunctionType>
416 auto dispatchOnBlocksize(FunctionType function) const
417 {
418 return m_matrix->dispatchOnBlocksize(function);
419 }
420
421private:
422 std::unique_ptr<matrix_type> m_matrix{};
423};
424
425} // namespace Opm::gpuistl
426
427#endif
The GpuSparseMatrixGeneric class uses cuSPARSE Generic API for sparse matrix operations.
Definition: GpuSparseMatrixGeneric.hpp:48
The GpuSparseMatrix class simple wrapper class for a CuSparse matrix.
Definition: GpuSparseMatrix.hpp:61
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
The GpuSparseMatrixWrapper Checks CUDA/HIP version and dispatches a version either using the old or t...
Definition: GpuSparseMatrixWrapper.hpp:62
GpuSparseMatrixWrapper(const GpuSparseMatrixWrapper &other)
Definition: GpuSparseMatrixWrapper.hpp:150
GpuVector< int > & getRowIndices()
getRowIndices returns the row indices used to represent the BSR structure.
Definition: GpuSparseMatrixWrapper.hpp:271
void updateNonzeroValues(const GpuSparseMatrixWrapper< T, OtherForceLegacy > &matrix)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
Definition: GpuSparseMatrixWrapper.hpp:389
virtual void mv(const GpuVector< T > &x, GpuVector< T > &y) const
mv performs matrix vector multiply y = Ax
Definition: GpuSparseMatrixWrapper.hpp:340
static constexpr int max_block_size
Maximum block size supported by this implementation.
Definition: GpuSparseMatrixWrapper.hpp:103
GpuVector< T > & getNonZeroValues()
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition: GpuSparseMatrixWrapper.hpp:251
matrix_type * operator->()
Definition: GpuSparseMatrixWrapper.hpp:82
const GpuVector< T > & getNonZeroValues() const
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition: GpuSparseMatrixWrapper.hpp:261
const matrix_type * operator->() const
Definition: GpuSparseMatrixWrapper.hpp:88
std::size_t N() const
N returns the number of rows (which is equal to the number of columns)
Definition: GpuSparseMatrixWrapper.hpp:232
static GpuSparseMatrixWrapper< T, ForceLegacy > 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:183
void setNonUnitDiagonal()
setNonUnitDiagonal sets the CuSparse flag that this has non-unit diagional.
Definition: GpuSparseMatrixWrapper.hpp:224
GpuSparseMatrixWrapper(const GpuVector< int > &rowIndices, const GpuVector< int > &columnIndices, std::size_t blockSize)
Definition: GpuSparseMatrixWrapper.hpp:143
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:378
virtual ~GpuSparseMatrixWrapper()=default
detail::GpuSparseMatrixDescription & getDescription()
getDescription the cusparse matrix description.
Definition: GpuSparseMatrixWrapper.hpp:330
void setUpperTriangular()
Definition: GpuSparseMatrixWrapper.hpp:194
std::size_t dim() const
dim returns the dimension of the vector space on which this matrix acts
Definition: GpuSparseMatrixWrapper.hpp:312
T field_type
Definition: GpuSparseMatrixWrapper.hpp:64
std::size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition: GpuSparseMatrixWrapper.hpp:241
GpuVector< int > & getColumnIndices()
getColumnIndices returns the column indices used to represent the BSR structure.
Definition: GpuSparseMatrixWrapper.hpp:291
void setUnitDiagonal()
setUnitDiagonal sets the CuSparse flag that this has unit diagional.
Definition: GpuSparseMatrixWrapper.hpp:214
GpuSparseMatrixWrapper(const T *nonZeroElements, const int *rowIndices, const int *columnIndices, std::size_t numberOfNonzeroBlocks, std::size_t blockSize, std::size_t numberOfRows)
Definition: GpuSparseMatrixWrapper.hpp:118
GpuSparseMatrixWrapper & operator=(const GpuSparseMatrixWrapper &)=delete
void setLowerTriangular()
setLowerTriangular sets the CuSparse flag that this is an lower diagonal (with non-unit diagonal) mat...
Definition: GpuSparseMatrixWrapper.hpp:204
virtual void usmv(T alpha, const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=alpha * Ax + y
Definition: GpuSparseMatrixWrapper.hpp:362
GpuSparseMatrix< T > matrix_type
Definition: GpuSparseMatrixWrapper.hpp:74
const GpuVector< int > & getColumnIndices() const
getColumnIndices returns the column indices used to represent the BSR structure.
Definition: GpuSparseMatrixWrapper.hpp:301
const matrix_type & get() const
Definition: GpuSparseMatrixWrapper.hpp:166
std::size_t blockSize() const
blockSize size of the blocks
Definition: GpuSparseMatrixWrapper.hpp:320
const GpuVector< int > & getRowIndices() const
getRowIndices returns the row indices used to represent the BSR structure.
Definition: GpuSparseMatrixWrapper.hpp:281
auto dispatchOnBlocksize(FunctionType function) const
Dispatches a function based on the block size of the matrix.
Definition: GpuSparseMatrixWrapper.hpp:416
virtual void umv(const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=Ax+y
Definition: GpuSparseMatrixWrapper.hpp:350
The CuSparseResource class wraps a CuSparse resource in a proper RAII pattern.
Definition: CuSparseResource.hpp:55
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition: AmgxInterface.hpp:38