19#ifndef OPM_GPUSPARSEMATRIXGENERIC_HPP
20#define OPM_GPUSPARSEMATRIXGENERIC_HPP
22#include <opm/common/ErrorMacros.hpp>
23#include <opm/simulators/linalg/gpuistl/GpuBuffer.hpp>
24#include <opm/simulators/linalg/gpuistl/GpuVector.hpp>
25#include <opm/simulators/linalg/gpuistl/detail/CuSparseHandle.hpp>
26#include <opm/simulators/linalg/gpuistl/detail/cusparse_safe_call.hpp>
27#include <opm/simulators/linalg/gpuistl/detail/gpusparse_matrix_utilities.hpp>
28#include <opm/simulators/linalg/gpuistl/detail/safe_conversion.hpp>
75 const int* rowIndices,
76 const int* columnIndices,
77 size_t numberOfNonzeroBlocks,
113 template <
class MatrixType>
145 return m_nonZeroElements;
155 return m_nonZeroElements;
185 return m_columnIndices;
195 return m_columnIndices;
260 template <
class MatrixType>
261 void updateNonzeroValues(
const MatrixType& matrix,
bool copyNonZeroElementsDirectly =
false);
290 template<
class FunctionType>
293 return dispatchOnBlocksizeImpl<max_block_size>(function);
305 const int m_numberOfNonzeroBlocks;
306 const int m_numberOfRows;
307 const int m_blockSize;
314 mutable GpuBuffer<std::byte> m_buffer;
320 void initializeMatrixDescriptor();
322 template <
class VectorType>
323 void assertSameSize(
const VectorType& vector)
const;
326 constexpr cudaDataType getDataType()
const
328 if constexpr (std::is_same_v<T, float>) {
330 }
else if constexpr (std::is_same_v<T, double>) {
333 static_assert(std::is_same_v<T, float> || std::is_same_v<T, double>,
"Only float and double are supported");
338 template<
int blockSizeCompileTime,
class FunctionType>
339 auto dispatchOnBlocksizeImpl(FunctionType function)
const
341 if (blockSizeCompileTime == m_blockSize) {
342 return function(std::integral_constant<int, blockSizeCompileTime>());
345 if constexpr (blockSizeCompileTime > 1) {
346 return dispatchOnBlocksizeImpl<blockSizeCompileTime - 1>(function);
348 OPM_THROW(std::runtime_error, fmt::format(
"Unsupported block size: {}", m_blockSize));
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:204
const GpuVector< int > & getColumnIndices() const
getColumnIndices returns the column indices used to represent the BSR structure.
Definition GpuSparseMatrixGeneric.hpp:193
size_t blockSize() const
blockSize size of the blocks
Definition GpuSparseMatrixGeneric.hpp:217
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 GpuSparseMatrixGeneric.cpp:231
GpuVector< int > & getColumnIndices()
getColumnIndices returns the column indices used to represent the BSR structure.
Definition GpuSparseMatrixGeneric.hpp:183
const GpuVector< T > & getNonZeroValues() const
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition GpuSparseMatrixGeneric.hpp:153
GpuVector< T > & getNonZeroValues()
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition GpuSparseMatrixGeneric.hpp:143
GpuVector< int > & getRowIndices()
getRowIndices returns the row indices used to represent the BSR structure.
Definition GpuSparseMatrixGeneric.hpp:163
size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition GpuSparseMatrixGeneric.hpp:128
virtual void umv(const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=Ax+y
Definition GpuSparseMatrixGeneric.cpp:299
virtual void mv(const GpuVector< T > &x, GpuVector< T > &y) const
mv performs matrix vector multiply y = Ax
Definition GpuSparseMatrixGeneric.cpp:290
GpuSparseMatrixGeneric(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 GpuSparseMatrixGeneric.cpp:42
void preprocessSpMV()
Preprocess SpMV operation to optimize for sparsity pattern.
Definition GpuSparseMatrixGeneric.cpp:140
static constexpr int max_block_size
Maximum block size supported by this implementation.
Definition GpuSparseMatrixGeneric.hpp:63
size_t N() const
N returns the number of rows (which is equal to the number of columns).
Definition GpuSparseMatrixGeneric.hpp:119
auto dispatchOnBlocksize(FunctionType function) const
Dispatches a function based on the block size of the matrix.
Definition GpuSparseMatrixGeneric.hpp:291
virtual void usmv(T alpha, const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=alpha * Ax + y
Definition GpuSparseMatrixGeneric.cpp:308
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
Definition GpuSparseMatrixGeneric.cpp:202
const GpuVector< int > & getRowIndices() const
getRowIndices returns the row indices used to represent the BSR structure.
Definition GpuSparseMatrixGeneric.hpp:173
Definition gpu_type_detection.hpp:30
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