19#ifndef OPM_GPUSPARSEMATRIX_HPP
20#define OPM_GPUSPARSEMATRIX_HPP
22#include <opm/common/ErrorMacros.hpp>
23#include <opm/simulators/linalg/gpuistl/detail/CuMatrixDescription.hpp>
24#include <opm/simulators/linalg/gpuistl/detail/CuSparseHandle.hpp>
25#include <opm/simulators/linalg/gpuistl/detail/safe_conversion.hpp>
26#include <opm/simulators/linalg/gpuistl/GpuVector.hpp>
27#include <opm/simulators/linalg/gpuistl/GpuSparseMatrixGeneric.hpp>
89 const int* rowIndices,
90 const int* columnIndices,
91 std::size_t numberOfNonzeroBlocks,
93 std::size_t numberOfRows);
126 template <
class MatrixType>
153 std::size_t
N()
const
184 if (m_genericMatrixForBlockSize1) {
185 return m_genericMatrixForBlockSize1->getNonZeroValues();
187 return m_nonZeroElements;
197 if (m_genericMatrixForBlockSize1) {
198 return m_genericMatrixForBlockSize1->getNonZeroValues();
200 return m_nonZeroElements;
210 if (m_genericMatrixForBlockSize1) {
211 return m_genericMatrixForBlockSize1->getRowIndices();
223 if (m_genericMatrixForBlockSize1) {
224 return m_genericMatrixForBlockSize1->getRowIndices();
236 if (m_genericMatrixForBlockSize1) {
237 return m_genericMatrixForBlockSize1->getColumnIndices();
239 return m_columnIndices;
249 if (m_genericMatrixForBlockSize1) {
250 return m_genericMatrixForBlockSize1->getColumnIndices();
252 return m_columnIndices;
291 return *m_matrixDescription;
327 template <
class MatrixType>
328 void updateNonzeroValues(
const MatrixType& matrix,
bool copyNonZeroElementsDirectly =
false);
369 template<
class FunctionType>
372 return dispatchOnBlocksizeImpl<max_block_size>(function);
384 const int m_numberOfNonzeroBlocks;
385 const int m_numberOfRows;
386 const int m_blockSize;
392 std::unique_ptr<GpuSparseMatrixGeneric<T>> m_genericMatrixForBlockSize1;
394 template <
class VectorType>
395 void assertSameSize(
const VectorType& vector)
const;
397 template<
int blockSizeCompileTime,
class FunctionType>
398 auto dispatchOnBlocksizeImpl(FunctionType function)
const
400 if (blockSizeCompileTime == m_blockSize) {
401 return function(std::integral_constant<int, blockSizeCompileTime>());
404 if constexpr (blockSizeCompileTime > 1) {
405 return dispatchOnBlocksizeImpl<blockSizeCompileTime - 1>(function);
407 OPM_THROW(std::runtime_error, fmt::format(fmt::runtime(
"Unsupported block size: {}"), m_blockSize));
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
const GpuVector< int > & getRowIndices() const
getRowIndices returns the row indices used to represent the BSR structure.
Definition GpuSparseMatrix.hpp:221
void setToZero()
setToZero resets the matrix to zero values.
Definition GpuSparseMatrix.cpp:207
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
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:116
void setNonUnitDiagonal()
setNonUnitDiagonal sets the CuSparse flag that this has non-unit diagional.
Definition GpuSparseMatrix.cpp:241
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 GpuSparseMatrix.cpp:149
virtual void umv(const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=Ax+y
Definition GpuSparseMatrix.cpp:283
GpuSparseMatrix(const GpuVector< int > &rowIndices, const GpuVector< int > &columnIndices, std::size_t blockSize)
Create a sparse matrix by copying the sparsity structure of another matrix, not filling in the values...
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
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 setUpperTriangular()
setUpperTriangular sets the CuSparse flag that this is an upper diagonal (with unit diagonal) matrix.
Definition GpuSparseMatrix.cpp:220
void setLowerTriangular()
setLowerTriangular sets the CuSparse flag that this is an lower diagonal (with non-unit diagonal) mat...
Definition GpuSparseMatrix.cpp:227
static constexpr int max_block_size
Maximum block size supported by this implementation.
Definition GpuSparseMatrix.hpp:73
virtual void mv(const GpuVector< T > &x, GpuVector< T > &y) const
mv performs matrix vector multiply y = Ax
Definition GpuSparseMatrix.cpp:248
virtual void usmv(T alpha, const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=alpha * Ax + y
Definition GpuSparseMatrix.cpp:319
auto dispatchOnBlocksize(FunctionType function) const
Dispatches a function based on the block size of the matrix.
Definition GpuSparseMatrix.hpp:370
void setUnitDiagonal()
setUnitDiagonal sets the CuSparse flag that this has unit diagional.
Definition GpuSparseMatrix.cpp:234
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)
Create the sparse matrix specified by the raw data.
Definition gpu_type_detection.hpp:30
The CuSparseHandle class provides a singleton for the simulator universal cuSparseHandle.
Definition CuSparseHandle.hpp:41
CuSparseResource< cusparseMatDescr_t > GpuSparseMatrixDescription
GpuSparseMatrixDescription holder.
Definition CuMatrixDescription.hpp:30
std::shared_ptr< CuSparseResource< cusparseMatDescr_t > > GpuSparseMatrixDescriptionPtr
Pointer to GpuSparseMatrixDescription holder.
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:97
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition AmgxInterface.hpp:38