GpuSparseMatrix.hpp
Go to the documentation of this file.
1/*
2 Copyright 2022-2023 SINTEF AS
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_GPUSPARSEMATRIX_HPP
20#define OPM_GPUSPARSEMATRIX_HPP
21
22#include <opm/common/ErrorMacros.hpp>
28
29#include <cstddef>
30#include <cusparse.h>
31#include <memory>
32#include <stdexcept>
33#include <type_traits>
34
35namespace Opm::gpuistl
36{
37
58template <typename T>
60
61{
62public:
63 using field_type = T;
64
73 static constexpr int max_block_size = 6;
74
88 GpuSparseMatrix(const T* nonZeroElements,
89 const int* rowIndices,
90 const int* columnIndices,
91 std::size_t numberOfNonzeroBlocks,
92 std::size_t blockSize,
93 std::size_t numberOfRows);
94
106 const GpuVector<int>& columnIndices,
107 std::size_t blockSize);
108
110
111 // We want to have this as non-mutable as possible, that is we do not want
112 // to deal with changing matrix sizes and sparsity patterns.
114
116
126 template <class MatrixType>
127 static GpuSparseMatrix<T> fromMatrix(const MatrixType& matrix, bool copyNonZeroElementsDirectly = false);
128
133
138
143
144
149
153 std::size_t N() const
154 {
155 // Technically this safe conversion is not needed since we enforce these to be
156 // non-negative in the constructor, but keeping them for added sanity for now.
157 //
158 // We don't believe this will yield any performance penality (it's used too far away from the inner loop),
159 // but should that be false, they can be removed.
160 return detail::to_size_t(m_numberOfRows);
161 }
162
167 std::size_t nonzeroes() const
168 {
169 // Technically this safe conversion is not needed since we enforce these to be
170 // non-negative in the constructor, but keeping them for added sanity for now.
171 //
172 // We don't believe this will yield any performance penality (it's used too far away from the inner loop),
173 // but should that be false, they can be removed.
174 return detail::to_size_t(m_numberOfNonzeroBlocks);
175 }
176
183 {
184 if (m_genericMatrixForBlockSize1) {
185 return m_genericMatrixForBlockSize1->getNonZeroValues();
186 }
187 return m_nonZeroElements;
188 }
189
196 {
197 if (m_genericMatrixForBlockSize1) {
198 return m_genericMatrixForBlockSize1->getNonZeroValues();
199 }
200 return m_nonZeroElements;
201 }
202
209 {
210 if (m_genericMatrixForBlockSize1) {
211 return m_genericMatrixForBlockSize1->getRowIndices();
212 }
213 return m_rowIndices;
214 }
215
222 {
223 if (m_genericMatrixForBlockSize1) {
224 return m_genericMatrixForBlockSize1->getRowIndices();
225 }
226 return m_rowIndices;
227 }
228
235 {
236 if (m_genericMatrixForBlockSize1) {
237 return m_genericMatrixForBlockSize1->getColumnIndices();
238 }
239 return m_columnIndices;
240 }
241
248 {
249 if (m_genericMatrixForBlockSize1) {
250 return m_genericMatrixForBlockSize1->getColumnIndices();
251 }
252 return m_columnIndices;
253 }
254
261 std::size_t dim() const
262 {
263 // Technically this safe conversion is not needed since we enforce these to be
264 // non-negative in the constructor, but keeping them for added sanity for now.
265 //
266 // We don't believe this will yield any performance penality (it's used too far away from the inner loop),
267 // but should that be false, they can be removed.
268 return detail::to_size_t(m_blockSize) * detail::to_size_t(m_numberOfRows);
269 }
270
274 std::size_t blockSize() const
275 {
276 // Technically this safe conversion is not needed since we enforce these to be
277 // non-negative in the constructor, but keeping them for added sanity for now.
278 //
279 // We don't believe this will yield any performance penality (it's used too far away from the inner loop),
280 // but should that be false, they can be removed.
281 return detail::to_size_t(m_blockSize);
282 }
283
290 {
291 return *m_matrixDescription;
292 }
293
299 virtual void mv(const GpuVector<T>& x, GpuVector<T>& y) const;
300
306 virtual void umv(const GpuVector<T>& x, GpuVector<T>& y) const;
307
308
315 virtual void usmv(T alpha, const GpuVector<T>& x, GpuVector<T>& y) const;
316
327 template <class MatrixType>
328 void updateNonzeroValues(const MatrixType& matrix, bool copyNonZeroElementsDirectly = false);
329
336
343
364 template<class FunctionType>
365 auto dispatchOnBlocksize(FunctionType function) const
366 {
367 return dispatchOnBlocksizeImpl<max_block_size>(function);
368 }
369
370private:
371 GpuVector<T> m_nonZeroElements;
372 GpuVector<int> m_columnIndices;
373 GpuVector<int> m_rowIndices;
374
375 // Notice that we store these three as int to make sure we are cusparse compatible.
376 //
377 // This gives the added benefit of checking the size constraints at construction of the matrix
378 // rather than in some call to cusparse.
379 const int m_numberOfNonzeroBlocks;
380 const int m_numberOfRows;
381 const int m_blockSize;
382
383 detail::GpuSparseMatrixDescriptionPtr m_matrixDescription;
384 detail::CuSparseHandle& m_cusparseHandle;
385
386 // For blockSize == 1, we use the generic API
387 std::unique_ptr<GpuSparseMatrixGeneric<T>> m_genericMatrixForBlockSize1;
388
389 template <class VectorType>
390 void assertSameSize(const VectorType& vector) const;
391
392 template<int blockSizeCompileTime, class FunctionType>
393 auto dispatchOnBlocksizeImpl(FunctionType function) const
394 {
395 if (blockSizeCompileTime == m_blockSize) {
396 return function(std::integral_constant<int, blockSizeCompileTime>());
397 }
398
399 if constexpr (blockSizeCompileTime > 1) {
400 return dispatchOnBlocksizeImpl<blockSizeCompileTime - 1>(function);
401 } else {
402 OPM_THROW(std::runtime_error, fmt::format("Unsupported block size: {}", m_blockSize));
403 }
404 }
405};
406
407} // namespace Opm::gpuistl
408
409#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
GpuSparseMatrix & operator=(const GpuSparseMatrix &)=delete
virtual void mv(const GpuVector< T > &x, GpuVector< T > &y) const
mv performs matrix vector multiply y = Ax
const GpuVector< int > & getRowIndices() const
getRowIndices returns the row indices used to represent the BSR structure.
Definition: GpuSparseMatrix.hpp:221
std::size_t blockSize() const
blockSize size of the blocks
Definition: GpuSparseMatrix.hpp:274
std::size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition: GpuSparseMatrix.hpp:167
const GpuVector< int > & getColumnIndices() const
getColumnIndices returns the column indices used to represent the BSR structure.
Definition: GpuSparseMatrix.hpp:247
void setNonUnitDiagonal()
setNonUnitDiagonal sets the CuSparse flag that this has non-unit diagional.
T field_type
Definition: GpuSparseMatrix.hpp:63
void updateNonzeroValues(const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
void updateNonzeroValues(const GpuSparseMatrix< T > &matrix)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
GpuSparseMatrix(const GpuVector< int > &rowIndices, const GpuVector< int > &columnIndices, std::size_t blockSize)
std::size_t dim() const
dim returns the dimension of the vector space on which this matrix acts
Definition: GpuSparseMatrix.hpp:261
detail::GpuSparseMatrixDescription & getDescription()
getDescription the cusparse matrix description.
Definition: GpuSparseMatrix.hpp:289
GpuSparseMatrix(const GpuSparseMatrix &)
GpuVector< int > & getRowIndices()
getRowIndices returns the row indices used to represent the BSR structure.
Definition: GpuSparseMatrix.hpp:208
const GpuVector< T > & getNonZeroValues() const
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition: GpuSparseMatrix.hpp:195
void updateNonzeroValues(const GpuSparseMatrixGeneric< T > &matrix)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
void setUpperTriangular()
setUpperTriangular sets the CuSparse flag that this is an upper diagonal (with unit diagonal) matrix.
virtual void umv(const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=Ax+y
void setLowerTriangular()
setLowerTriangular sets the CuSparse flag that this is an lower diagonal (with non-unit diagonal) mat...
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
static constexpr int max_block_size
Maximum block size supported by this implementation.
Definition: GpuSparseMatrix.hpp:73
auto dispatchOnBlocksize(FunctionType function) const
Dispatches a function based on the block size of the matrix.
Definition: GpuSparseMatrix.hpp:365
void setUnitDiagonal()
setUnitDiagonal sets the CuSparse flag that this has unit diagional.
virtual void usmv(T alpha, const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=alpha * Ax + y
GpuVector< int > & getColumnIndices()
getColumnIndices returns the column indices used to represent the BSR structure.
Definition: GpuSparseMatrix.hpp:234
std::size_t N() const
N returns the number of rows (which is equal to the number of columns)
Definition: GpuSparseMatrix.hpp:153
GpuVector< T > & getNonZeroValues()
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition: GpuSparseMatrix.hpp:182
GpuSparseMatrix(const T *nonZeroElements, const int *rowIndices, const int *columnIndices, std::size_t numberOfNonzeroBlocks, std::size_t blockSize, std::size_t numberOfRows)
The CuSparseHandle class provides a singleton for the simulator universal cuSparseHandle.
Definition: CuSparseHandle.hpp:41
The CuSparseResource class wraps a CuSparse resource in a proper RAII pattern.
Definition: CuSparseResource.hpp:55
std::shared_ptr< CuSparseResource< cusparseMatDescr_t > > GpuSparseMatrixDescriptionPtr
Definition: CuMatrixDescription.hpp:35
__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
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition: AmgxInterface.hpp:38