Opm::cuistl::CuSparseMatrix< T > Class Template Reference

The CuSparseMatrix class simple wrapper class for a CuSparse matrix. More...

#include <CuSparseMatrix.hpp>

Inheritance diagram for Opm::cuistl::CuSparseMatrix< T >:
Inheritance graph

Public Member Functions

 CuSparseMatrix (const T *nonZeroElements, const int *rowIndices, const int *columnIndices, size_t numberOfNonzeroBlocks, size_t blockSize, size_t numberOfRows)
 
 CuSparseMatrix (const CuSparseMatrix &)=delete
 
CuSparseMatrixoperator= (const CuSparseMatrix &)=delete
 
virtual ~CuSparseMatrix ()
 
void setUpperTriangular ()
 setUpperTriangular sets the CuSparse flag that this is an upper diagonal (with unit diagonal) matrix. More...
 
void setLowerTriangular ()
 setLowerTriangular sets the CuSparse flag that this is an lower diagonal (with non-unit diagonal) matrix. More...
 
void setUnitDiagonal ()
 setUnitDiagonal sets the CuSparse flag that this has unit diagional. More...
 
void setNonUnitDiagonal ()
 setNonUnitDiagonal sets the CuSparse flag that this has non-unit diagional. More...
 
size_t N () const
 N returns the number of rows (which is equal to the number of columns) More...
 
size_t nonzeroes () const
 nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blocks More...
 
CuVector< T > & getNonZeroValues ()
 getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block) More...
 
const CuVector< T > & getNonZeroValues () const
 getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block) More...
 
CuVector< int > & getRowIndices ()
 getRowIndices returns the row indices used to represent the BSR structure. More...
 
const CuVector< int > & getRowIndices () const
 getRowIndices returns the row indices used to represent the BSR structure. More...
 
CuVector< int > & getColumnIndices ()
 getColumnIndices returns the column indices used to represent the BSR structure. More...
 
const CuVector< int > & getColumnIndices () const
 getColumnIndices returns the column indices used to represent the BSR structure. More...
 
size_t dim () const
 dim returns the dimension of the vector space on which this matrix acts More...
 
size_t blockSize () const
 blockSize size of the blocks More...
 
detail::CuSparseMatrixDescriptiongetDescription ()
 getDescription the cusparse matrix description. More...
 
virtual void mv (const CuVector< T > &x, CuVector< T > &y) const
 mv performs matrix vector multiply y = Ax More...
 
virtual void umv (const CuVector< T > &x, CuVector< T > &y) const
 umv computes y=Ax+y More...
 
virtual void usmv (T alpha, const CuVector< T > &x, CuVector< T > &y) const
 umv computes y=alpha * Ax + y More...
 
template<class MatrixType >
void updateNonzeroValues (const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
 updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix More...
 

Static Public Member Functions

template<class MatrixType >
static CuSparseMatrix< T > fromMatrix (const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
 fromMatrix creates a new matrix with the same block size and values as the given matrix More...
 

Detailed Description

template<typename T>
class Opm::cuistl::CuSparseMatrix< T >

The CuSparseMatrix class simple wrapper class for a CuSparse matrix.

Note
we currently only support simple raw primitives for T (double and float). Block size is handled through the block size parameter
Template Parameters
Tthe type to store. Can be either float, double or int.
Note
we only support square matrices.
We only support Block Compressed Sparse Row Format (BSR) for now.

Constructor & Destructor Documentation

◆ CuSparseMatrix() [1/2]

template<typename T >
Opm::cuistl::CuSparseMatrix< T >::CuSparseMatrix ( 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.

Note
Prefer to use the constructor taking a const reference to a matrix instead.
Parameters
[in]nonZeroElementsthe non-zero values of the matrix
[in]rowIndicesthe row indices of the non-zero elements
[in]columnIndicesthe column indices of the non-zero elements
[in]numberOfNonzeroElementsnumber of nonzero elements
[in]blockSizesize of each block matrix (typically 3)
[in]numberOfRowsthe number of rows
Note
We assume numberOfNonzeroBlocks, blockSize and numberOfRows all are representable as int due to restrictions in the current version of cusparse. This might change in future versions.

◆ CuSparseMatrix() [2/2]

template<typename T >
Opm::cuistl::CuSparseMatrix< T >::CuSparseMatrix ( const CuSparseMatrix< T > &  )
delete

We don't want to be able to copy this for now (too much hassle in copying the cusparse resources)

◆ ~CuSparseMatrix()

template<typename T >
virtual Opm::cuistl::CuSparseMatrix< T >::~CuSparseMatrix ( )
virtual

Member Function Documentation

◆ blockSize()

template<typename T >
size_t Opm::cuistl::CuSparseMatrix< T >::blockSize ( ) const
inline

blockSize size of the blocks

References Opm::cuistl::detail::to_size_t().

◆ dim()

template<typename T >
size_t Opm::cuistl::CuSparseMatrix< T >::dim ( ) const
inline

dim returns the dimension of the vector space on which this matrix acts

This is equivalent to matrix.N() * matrix.blockSize()

Returns
matrix.N() * matrix.blockSize()

References Opm::cuistl::detail::to_size_t().

◆ fromMatrix()

template<typename T >
template<class MatrixType >
static CuSparseMatrix< T > Opm::cuistl::CuSparseMatrix< T >::fromMatrix ( const MatrixType &  matrix,
bool  copyNonZeroElementsDirectly = false 
)
static

fromMatrix creates a new matrix with the same block size and values as the given matrix

Parameters
matrixthe matrix to copy from
copyNonZeroElementsDirectlyif true will do a memcpy from matrix[0][0][0][0], otherwise will build up the non-zero elements by looping over the matrix. Note that setting this to true will yield a performance increase, but might not always yield correct results depending on how the matrix has been initialized. If unsure, leave it as false.
Template Parameters
MatrixTypeis assumed to be a Dune::BCRSMatrix compatible matrix.

◆ getColumnIndices() [1/2]

template<typename T >
CuVector< int > & Opm::cuistl::CuSparseMatrix< T >::getColumnIndices ( )
inline

getColumnIndices returns the column indices used to represent the BSR structure.

Returns
Read the CuSPARSE documentation on Block Compressed Sparse Row Format (BSR) for the exact ordering.

◆ getColumnIndices() [2/2]

template<typename T >
const CuVector< int > & Opm::cuistl::CuSparseMatrix< T >::getColumnIndices ( ) const
inline

getColumnIndices returns the column indices used to represent the BSR structure.

Returns
Read the CuSPARSE documentation on Block Compressed Sparse Row Format (BSR) for the exact ordering.

◆ getDescription()

template<typename T >
detail::CuSparseMatrixDescription & Opm::cuistl::CuSparseMatrix< T >::getDescription ( )
inline

getDescription the cusparse matrix description.

This description is needed for most calls to the CuSparse library

◆ getNonZeroValues() [1/2]

template<typename T >
CuVector< T > & Opm::cuistl::CuSparseMatrix< T >::getNonZeroValues ( )
inline

getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)

Note
Read the CuSPARSE documentation on Block Compressed Sparse Row Format (BSR) for the exact ordering.

◆ getNonZeroValues() [2/2]

template<typename T >
const CuVector< T > & Opm::cuistl::CuSparseMatrix< T >::getNonZeroValues ( ) const
inline

getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)

Note
Read the CuSPARSE documentation on Block Compressed Sparse Row Format (BSR) for the exact ordering.

◆ getRowIndices() [1/2]

template<typename T >
CuVector< int > & Opm::cuistl::CuSparseMatrix< T >::getRowIndices ( )
inline

getRowIndices returns the row indices used to represent the BSR structure.

Note
Read the CuSPARSE documentation on Block Compressed Sparse Row Format (BSR) for the exact ordering.

◆ getRowIndices() [2/2]

template<typename T >
const CuVector< int > & Opm::cuistl::CuSparseMatrix< T >::getRowIndices ( ) const
inline

getRowIndices returns the row indices used to represent the BSR structure.

Note
Read the CuSPARSE documentation on Block Compressed Sparse Row Format (BSR) for the exact ordering.

◆ mv()

template<typename T >
virtual void Opm::cuistl::CuSparseMatrix< T >::mv ( const CuVector< T > &  x,
CuVector< T > &  y 
) const
virtual

mv performs matrix vector multiply y = Ax

Parameters
[in]xthe vector to multiply the matrix with
[out]ythe output vector
Note
Due to limitations of CuSparse, this is only supported for block sizes greater than 1.

◆ N()

template<typename T >
size_t Opm::cuistl::CuSparseMatrix< T >::N ( ) const
inline

N returns the number of rows (which is equal to the number of columns)

References Opm::cuistl::detail::to_size_t().

◆ nonzeroes()

template<typename T >
size_t Opm::cuistl::CuSparseMatrix< T >::nonzeroes ( ) const
inline

nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blocks

Returns
number of non zero blocks.

References Opm::cuistl::detail::to_size_t().

◆ operator=()

template<typename T >
CuSparseMatrix & Opm::cuistl::CuSparseMatrix< T >::operator= ( const CuSparseMatrix< T > &  )
delete

We don't want to be able to copy this for now (too much hassle in copying the cusparse resources)

◆ setLowerTriangular()

template<typename T >
void Opm::cuistl::CuSparseMatrix< T >::setLowerTriangular ( )

setLowerTriangular sets the CuSparse flag that this is an lower diagonal (with non-unit diagonal) matrix.

◆ setNonUnitDiagonal()

template<typename T >
void Opm::cuistl::CuSparseMatrix< T >::setNonUnitDiagonal ( )

setNonUnitDiagonal sets the CuSparse flag that this has non-unit diagional.

◆ setUnitDiagonal()

template<typename T >
void Opm::cuistl::CuSparseMatrix< T >::setUnitDiagonal ( )

setUnitDiagonal sets the CuSparse flag that this has unit diagional.

◆ setUpperTriangular()

template<typename T >
void Opm::cuistl::CuSparseMatrix< T >::setUpperTriangular ( )

setUpperTriangular sets the CuSparse flag that this is an upper diagonal (with unit diagonal) matrix.

◆ umv()

template<typename T >
virtual void Opm::cuistl::CuSparseMatrix< T >::umv ( const CuVector< T > &  x,
CuVector< T > &  y 
) const
virtual

umv computes y=Ax+y

Parameters
[in]xthe vector to multiply with A
[in,out]ythe vector to add and store the output in
Note
Due to limitations of CuSparse, this is only supported for block sizes greater than 1.

◆ updateNonzeroValues()

template<typename T >
template<class MatrixType >
void Opm::cuistl::CuSparseMatrix< T >::updateNonzeroValues ( const MatrixType &  matrix,
bool  copyNonZeroElementsDirectly = false 
)

updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix

Parameters
matrixthe matrix to extract the non-zero values from
copyNonZeroElementsDirectlyif true will do a memcpy from matrix[0][0][0][0], otherwise will build up the non-zero elements by looping over the matrix. Note that setting this to true will yield a performance increase, but might not always yield correct results depending on how the matrix matrix has been initialized. If unsure, leave it as false.
Note
This assumes the given matrix has the same sparsity pattern.
Template Parameters
MatrixTypeis assumed to be a Dune::BCRSMatrix compatible matrix.

Referenced by Opm::cuistl::SolverAdapter< Operator, UnderlyingSolver, X >::apply().

◆ usmv()

template<typename T >
virtual void Opm::cuistl::CuSparseMatrix< T >::usmv ( alpha,
const CuVector< T > &  x,
CuVector< T > &  y 
) const
virtual

umv computes y=alpha * Ax + y

Parameters
[in]xthe vector to multiply with A
[in,out]ythe vector to add and store the output in
Note
Due to limitations of CuSparse, this is only supported for block sizes greater than 1.

The documentation for this class was generated from the following file: