19#ifndef OPM_GPUSPARSEMATRIXGENERIC_HPP
20#define OPM_GPUSPARSEMATRIXGENERIC_HPP
22#include <opm/common/ErrorMacros.hpp>
23#include <opm/simulators/linalg/gpuistl/detail/CuSparseHandle.hpp>
24#include <opm/simulators/linalg/gpuistl/detail/cusparse_safe_call.hpp>
25#include <opm/simulators/linalg/gpuistl/detail/gpusparse_matrix_utilities.hpp>
26#include <opm/simulators/linalg/gpuistl/detail/safe_conversion.hpp>
27#include <opm/simulators/linalg/gpuistl/GpuBuffer.hpp>
28#include <opm/simulators/linalg/gpuistl/GpuVector.hpp>
72 const int* rowIndices,
73 const int* columnIndices,
74 std::size_t numberOfNonzeroBlocks,
76 std::size_t numberOfRows);
110 template <
class MatrixType>
116 std::size_t
N()
const
142 return m_nonZeroElements;
152 return m_nonZeroElements;
182 return m_columnIndices;
192 return m_columnIndices;
257 template <
class MatrixType>
258 void updateNonzeroValues(
const MatrixType& matrix,
bool copyNonZeroElementsDirectly =
false);
292 template<
class FunctionType>
295 return dispatchOnBlocksizeImpl<max_block_size>(function);
307 const int m_numberOfNonzeroBlocks;
308 const int m_numberOfRows;
309 const int m_blockSize;
316 mutable GpuBuffer<std::byte> m_buffer;
322 void initializeMatrixDescriptor();
324 template <
class VectorType>
325 void assertSameSize(
const VectorType& vector)
const;
328 constexpr cudaDataType getDataType()
const
330 if constexpr (std::is_same_v<T, float>) {
332 }
else if constexpr (std::is_same_v<T, double>) {
335 static_assert(std::is_same_v<T, float> || std::is_same_v<T, double>,
"Only float and double are supported");
340 template<
int blockSizeCompileTime,
class FunctionType>
341 auto dispatchOnBlocksizeImpl(FunctionType function)
const
343 if (blockSizeCompileTime == m_blockSize) {
344 return function(std::integral_constant<int, blockSizeCompileTime>());
347 if constexpr (blockSizeCompileTime > 1) {
348 return dispatchOnBlocksizeImpl<blockSizeCompileTime - 1>(function);
350 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
std::size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition GpuSparseMatrixGeneric.hpp:125
std::size_t blockSize() const
blockSize size of the blocks
Definition GpuSparseMatrixGeneric.hpp:214
GpuSparseMatrixGeneric(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.
const GpuVector< int > & getColumnIndices() const
getColumnIndices returns the column indices used to represent the BSR structure.
Definition GpuSparseMatrixGeneric.hpp:190
void setToZero()
setToZero resets the matrix to zero values.
Definition GpuSparseMatrixGeneric.cpp:258
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:180
const GpuVector< T > & getNonZeroValues() const
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition GpuSparseMatrixGeneric.hpp:150
GpuVector< T > & getNonZeroValues()
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
Definition GpuSparseMatrixGeneric.hpp:140
GpuVector< int > & getRowIndices()
getRowIndices returns the row indices used to represent the BSR structure.
Definition GpuSparseMatrixGeneric.hpp:160
virtual void umv(const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=Ax+y
Definition GpuSparseMatrixGeneric.cpp:305
std::size_t dim() const
dim returns the dimension of the vector space on which this matrix acts
Definition GpuSparseMatrixGeneric.hpp:201
virtual void mv(const GpuVector< T > &x, GpuVector< T > &y) const
mv performs matrix vector multiply y = Ax
Definition GpuSparseMatrixGeneric.cpp:296
void preprocessSpMV()
Preprocess SpMV operation to optimize for sparsity pattern.
Definition GpuSparseMatrixGeneric.cpp:140
GpuSparseMatrixGeneric(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...
static constexpr int max_block_size
Maximum block size supported by this implementation.
Definition GpuSparseMatrixGeneric.hpp:60
std::size_t N() const
N returns the number of rows (which is equal to the number of columns).
Definition GpuSparseMatrixGeneric.hpp:116
auto dispatchOnBlocksize(FunctionType function) const
Dispatches a function based on the block size of the matrix.
Definition GpuSparseMatrixGeneric.hpp:293
virtual void usmv(T alpha, const GpuVector< T > &x, GpuVector< T > &y) const
umv computes y=alpha * Ax + y
Definition GpuSparseMatrixGeneric.cpp:314
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:170
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:97
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition AmgxInterface.hpp:38