gpusparse_matrix_utilities.hpp File Reference
#include <opm/common/ErrorMacros.hpp>
#include <opm/simulators/linalg/gpuistl/detail/cusparse_safe_call.hpp>
#include <opm/simulators/linalg/gpuistl/detail/safe_conversion.hpp>
#include <fmt/core.h>
#include <memory>
#include <vector>
#include <cusparse.h>
Include dependency graph for gpusparse_matrix_utilities.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Namespaces

namespace  Opm
 
namespace  Opm::gpuistl
 
namespace  Opm::gpuistl::detail
 

Macros

#define INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, realtype, blockdim)
 Macro for generating template instantiations for DUNE matrix operations. More...
 
#define INSTANTIATE_FOR_TYPE_AND_CLASS(CLASS_NAME, T)
 Macro for generating template instantiations for a range of block sizes. More...
 

Functions

auto Opm::gpuistl::detail::makeSafeVectorDescriptor ()
 Create RAII-managed cuSPARSE dense vector descriptor. More...
 
auto Opm::gpuistl::detail::makeSafeMatrixDescriptor ()
 Create RAII-managed cuSPARSE sparse matrix descriptor. More...
 
template<class T , class M >
std::vector< T > Opm::gpuistl::detail::extractNonzeroValues (const M &matrix)
 Extract non-zero values from a DUNE matrix in the order expected by cuSPARSE. More...
 
template<class M >
std::pair< std::vector< int >, std::vector< int > > Opm::gpuistl::detail::extractSparsityPattern (const M &matrix)
 Extract sparsity pattern from a DUNE matrix. More...
 
template<class T , class M >
T * Opm::gpuistl::detail::findFirstElementPointer (const M &matrix)
 Find pointer to first matrix element for direct memory access. More...
 
template<class M >
void Opm::gpuistl::detail::validateMatrixConversion (const M &matrix, const std::vector< int > &columnIndices, const std::vector< int > &rowIndices)
 Validate matrix conversion parameters. More...
 
void Opm::gpuistl::detail::validateMatrixCompatibility (size_t nonzeroes1, size_t blockSize1, size_t N1, size_t nonzeroes2, size_t blockSize2, size_t N2)
 Validate that two matrices have compatible dimensions. More...
 
void Opm::gpuistl::detail::validateVectorMatrixSizes (size_t vectorSize, size_t matrixBlockSize, size_t matrixColumns)
 Validate vector-matrix size compatibility for operations. More...
 

Macro Definition Documentation

◆ INSTANTIATE_FOR_TYPE_AND_CLASS

#define INSTANTIATE_FOR_TYPE_AND_CLASS (   CLASS_NAME,
 
)
Value:
INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 2); \
INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 3); \
INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 4); \
INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 5); \
INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 6)
#define INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, realtype, blockdim)
Macro for generating template instantiations for DUNE matrix operations.
Definition: gpusparse_matrix_utilities.hpp:231

Macro for generating template instantiations for a range of block sizes.

◆ INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS

#define INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS (   CLASS_NAME,
  realtype,
  blockdim 
)
Value:
template CLASS_NAME<realtype> CLASS_NAME<realtype>::fromMatrix( \
const Dune::BCRSMatrix<Dune::FieldMatrix<realtype, blockdim, blockdim>>&, bool); \
template CLASS_NAME<realtype> CLASS_NAME<realtype>::fromMatrix( \
const Dune::BCRSMatrix<Opm::MatrixBlock<realtype, blockdim, blockdim>>&, bool); \
template void CLASS_NAME<realtype>::updateNonzeroValues( \
const Dune::BCRSMatrix<Dune::FieldMatrix<realtype, blockdim, blockdim>>&, bool); \
template void CLASS_NAME<realtype>::updateNonzeroValues( \
const Dune::BCRSMatrix<Opm::MatrixBlock<realtype, blockdim, blockdim>>&, bool)
Definition: matrixblock.hh:227

Macro for generating template instantiations for DUNE matrix operations.

This macro generates the necessary explicit template instantiations for fromMatrix and updateNonzeroValues functions with DUNE matrix types.