GpuSparseMatrixGeneric.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_GPUSPARSEMATRIXGENERIC_HPP
20#define OPM_GPUSPARSEMATRIXGENERIC_HPP
21
22#include <opm/common/ErrorMacros.hpp>
29
30#include <cusparse.h>
31
32#include <cstddef>
33
34namespace Opm::gpuistl
35{
36
37
38
49template <typename T>
51{
52public:
53 using field_type = T;
54
64 GpuSparseMatrixGeneric(const T* nonZeroElements,
65 const int* rowIndices,
66 const int* columnIndices,
67 size_t numberOfNonzeroBlocks,
68 size_t blockSize,
69 size_t numberOfRows);
70
77 GpuSparseMatrixGeneric(const GpuVector<int>& rowIndices, const GpuVector<int>& columnIndices, size_t blockSize);
78
80
87
88 // We want to have this as non-mutable as possible, that is we do not want
89 // to deal with changing matrix sizes and sparsity patterns.
91
93
103 template <class MatrixType>
104 static GpuSparseMatrixGeneric<T> fromMatrix(const MatrixType& matrix, bool copyNonZeroElementsDirectly = false);
105
109 size_t N() const
110 {
111 return detail::to_size_t(m_numberOfRows);
112 }
113
118 size_t nonzeroes() const
119 {
120 // Technically this safe conversion is not needed since we enforce these to be
121 // non-negative in the constructor, but keeping them for added sanity for now.
122 //
123 // We don't believe this will yield any performance penality (it's used too far away from the inner loop),
124 // but should that be false, they can be removed.
125 return detail::to_size_t(m_numberOfNonzeroBlocks);
126 }
127
134 {
135 return m_nonZeroElements;
136 }
137
144 {
145 return m_nonZeroElements;
146 }
147
154 {
155 return m_rowIndices;
156 }
157
164 {
165 return m_rowIndices;
166 }
167
174 {
175 return m_columnIndices;
176 }
177
184 {
185 return m_columnIndices;
186 }
187
194 size_t dim() const
195 {
196 // Technically this safe conversion is not needed since we enforce these to be
197 // non-negative in the constructor, but keeping them for added sanity for now.
198 //
199 // We don't believe this will yield any performance penality (it's used too far away from the inner loop),
200 // but should that be false, they can be removed.
201 return detail::to_size_t(m_blockSize) * detail::to_size_t(m_numberOfRows);
202 }
203
207 size_t blockSize() const
208 {
209 // Technically this safe conversion is not needed since we enforce these to be
210 // non-negative in the constructor, but keeping them for added sanity for now.
211 //
212 // We don't believe this will yield any performance penality (it's used too far away from the inner loop),
213 // but should that be false, they can be removed.
214 return detail::to_size_t(m_blockSize);
215 }
216
222 virtual void mv(const GpuVector<T>& x, GpuVector<T>& y) const;
223
229 virtual void umv(const GpuVector<T>& x, GpuVector<T>& y) const;
230
231
238 virtual void usmv(T alpha, const GpuVector<T>& x, GpuVector<T>& y) const;
239
250 template <class MatrixType>
251 void updateNonzeroValues(const MatrixType& matrix, bool copyNonZeroElementsDirectly = false);
252
259
260private:
261 GpuVector<T> m_nonZeroElements;
262 GpuVector<int> m_columnIndices;
263 GpuVector<int> m_rowIndices;
264
265 // Notice that we store these three as int to make sure we are cusparse compatible.
266 //
267 // This gives the added benefit of checking the size constraints at construction of the matrix
268 // rather than in some call to cusparse.
269 const int m_numberOfNonzeroBlocks;
270 const int m_numberOfRows;
271 const int m_blockSize;
272
273 // Generic API descriptors
274 decltype(detail::makeSafeMatrixDescriptor()) m_matrixDescriptor;
275 detail::CuSparseHandle& m_cusparseHandle;
276
277 // Cached buffer for operations
278 mutable GpuBuffer<std::byte> m_buffer;
279
280 // Helper methods for SpMV operations
281 void spMV(T alpha, const GpuVector<T>& x, T beta, GpuVector<T>& y) const;
282
283 // Initialize matrix descriptor based on block size
284 void initializeMatrixDescriptor();
285
286 template <class VectorType>
287 void assertSameSize(const VectorType& vector) const;
288
289 // Helper to get cuSPARSE data type from C++ type
290 constexpr cudaDataType getDataType() const
291 {
292 if constexpr (std::is_same_v<T, float>) {
293 return CUDA_R_32F;
294 } else if constexpr (std::is_same_v<T, double>) {
295 return CUDA_R_64F;
296 } else {
297 static_assert(std::is_same_v<T, float> || std::is_same_v<T, double>, "Only float and double are supported");
298 return CUDA_R_32F; // Unreachable, but needed to compile
299 }
300 }
301};
302} // namespace Opm::gpuistl
303#endif
The GpuSparseMatrixGeneric class uses cuSPARSE Generic API for sparse matrix operations.
Definition: GpuSparseMatrixGeneric.hpp:51
size_t dim() const
dim returns the dimension of the vector space on which this matrix acts
Definition: GpuSparseMatrixGeneric.hpp:194
GpuSparseMatrixGeneric & operator=(const GpuSparseMatrixGeneric &)=delete
GpuSparseMatrixGeneric(const GpuVector< int > &rowIndices, const GpuVector< int > &columnIndices, size_t blockSize)
const GpuVector< int > & getColumnIndices() const
getColumnIndices returns the column indices used to represent the BSR structure.
Definition: GpuSparseMatrixGeneric.hpp:183
size_t blockSize() const
blockSize size of the blocks
Definition: GpuSparseMatrixGeneric.hpp:207
virtual void umv(const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=Ax+y
virtual void mv(const GpuVector< T > &x, GpuVector< T > &y) const
mv performs matrix vector multiply y = Ax
GpuSparseMatrixGeneric(const GpuSparseMatrixGeneric &)
void updateNonzeroValues(const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
GpuVector< int > & getColumnIndices()
getColumnIndices returns the column indices used to represent the BSR structure.
Definition: GpuSparseMatrixGeneric.hpp:173
const GpuVector< T > & getNonZeroValues() const
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition: GpuSparseMatrixGeneric.hpp:143
T field_type
Definition: GpuSparseMatrixGeneric.hpp:53
GpuVector< T > & getNonZeroValues()
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition: GpuSparseMatrixGeneric.hpp:133
void updateNonzeroValues(const GpuSparseMatrixGeneric< T > &matrix)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
GpuVector< int > & getRowIndices()
getRowIndices returns the row indices used to represent the BSR structure.
Definition: GpuSparseMatrixGeneric.hpp:153
size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition: GpuSparseMatrixGeneric.hpp:118
GpuSparseMatrixGeneric(const T *nonZeroElements, const int *rowIndices, const int *columnIndices, size_t numberOfNonzeroBlocks, size_t blockSize, size_t numberOfRows)
void preprocessSpMV()
Preprocess SpMV operation to optimize for sparsity pattern.
static GpuSparseMatrixGeneric< T > fromMatrix(const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
fromMatrix creates a new matrix with the same block size and values as the given matrix
size_t N() const
N returns the number of rows (which is equal to the number of columns)
Definition: GpuSparseMatrixGeneric.hpp:109
const GpuVector< int > & getRowIndices() const
getRowIndices returns the row indices used to represent the BSR structure.
Definition: GpuSparseMatrixGeneric.hpp:163
virtual void usmv(T alpha, const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=alpha * Ax + y
The CuSparseHandle class provides a singleton for the simulator universal cuSparseHandle.
Definition: CuSparseHandle.hpp:41
auto makeSafeMatrixDescriptor()
Create RAII-managed cuSPARSE sparse matrix descriptor.
Definition: gpusparse_matrix_utilities.hpp:56
__host__ __device__ std::size_t to_size_t(int i)
to_size_t converts a (on most relevant platforms) a 32 bit signed int to a 64 bits unsigned int
Definition: safe_conversion.hpp:86
Definition: AmgxInterface.hpp:38