|
opm-simulators
|
The GpuSparseMatrixGeneric class uses cuSPARSE Generic API for sparse matrix operations. More...
#include <GpuSparseMatrixGeneric.hpp>
Public Types | |
| using | field_type = T |
Public Member Functions | |
| 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. | |
| 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. | |
| GpuSparseMatrixGeneric (const GpuSparseMatrixGeneric &) | |
| void | preprocessSpMV () |
| Preprocess SpMV operation to optimize for sparsity pattern. | |
| GpuSparseMatrixGeneric & | operator= (const GpuSparseMatrixGeneric &)=delete |
| std::size_t | N () const |
| N returns the number of rows (which is equal to the number of columns). | |
| std::size_t | nonzeroes () const |
| nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blocks | |
| GpuVector< T > & | getNonZeroValues () |
| getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block) | |
| const GpuVector< T > & | getNonZeroValues () const |
| getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block) | |
| GpuVector< int > & | getRowIndices () |
| getRowIndices returns the row indices used to represent the BSR structure. | |
| const GpuVector< int > & | getRowIndices () const |
| getRowIndices returns the row indices used to represent the BSR structure. | |
| GpuVector< int > & | getColumnIndices () |
| getColumnIndices returns the column indices used to represent the BSR structure. | |
| const GpuVector< int > & | getColumnIndices () const |
| getColumnIndices returns the column indices used to represent the BSR structure. | |
| std::size_t | dim () const |
| dim returns the dimension of the vector space on which this matrix acts | |
| std::size_t | blockSize () const |
| blockSize size of the blocks | |
| virtual void | mv (const GpuVector< T > &x, GpuVector< T > &y) const |
| mv performs matrix vector multiply y = Ax | |
| virtual void | umv (const GpuVector< T > &x, GpuVector< T > &y) const |
| umv computes y=Ax+y | |
| virtual void | usmv (T alpha, const GpuVector< T > &x, GpuVector< T > &y) const |
| umv computes y=alpha * Ax + y | |
| 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 | |
| void | updateNonzeroValues (const GpuSparseMatrixGeneric< T > &matrix) |
| updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix | |
| void | setToZero () |
| setToZero resets the matrix to zero values. | |
| template<class FunctionType> | |
| auto | dispatchOnBlocksize (FunctionType function) const |
| Dispatches a function based on the block size of the matrix. | |
Static Public Member Functions | |
| template<class MatrixType> | |
| 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 | |
Static Public Attributes | |
| static constexpr int | max_block_size = 6 |
| Maximum block size supported by this implementation. | |
The GpuSparseMatrixGeneric class uses cuSPARSE Generic API for sparse matrix operations.
| T | the type to store. Can be either float or double. |
| Opm::gpuistl::GpuSparseMatrixGeneric< T >::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.
| [in] | nonZeroElements | the non-zero values of the matrix |
| [in] | rowIndices | the row indices of the non-zero elements |
| [in] | columnIndices | the column indices of the non-zero elements |
| [in] | numberOfNonzeroBlocks | number of nonzero elements |
| [in] | blockSize | size of each block matrix (typically 3) |
| [in] | numberOfRows | the number of rows |
| Opm::gpuistl::GpuSparseMatrixGeneric< T >::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.
| [in] | rowIndices | the row indices of the non-zero elements |
| [in] | columnIndices | the column indices of the non-zero elements |
| [in] | blockSize | size of each block matrix (typically 3) |
|
inline |
dim returns the dimension of the vector space on which this matrix acts
This is equivalent to matrix.N() * matrix.blockSize()
|
inline |
Dispatches a function based on the block size of the matrix.
This method allows executing different code paths depending on the block size of the matrix, up to the maximum block size specified by max_block_size.
Use this function if you need the block size to be known at compile time.
| FunctionType | Type of the function to be dispatched |
| function | The function to be executed based on the block size |
You can use this function as
|
static |
fromMatrix creates a new matrix with the same block size and values as the given matrix
| matrix | the matrix to copy from |
| copyNonZeroElementsDirectly | if 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. |
| MatrixType | is assumed to be a Dune::BCRSMatrix compatible matrix. |
|
inline |
getColumnIndices returns the column indices used to represent the BSR structure.
|
inline |
getColumnIndices returns the column indices used to represent the BSR structure.
|
inline |
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
|
inline |
getNonZeroValues returns the GPU vector containing the non-zero values (ordered by block)
|
inline |
getRowIndices returns the row indices used to represent the BSR structure.
|
inline |
getRowIndices returns the row indices used to represent the BSR structure.
|
virtual |
mv performs matrix vector multiply y = Ax
| [in] | x | the vector to multiply the matrix with |
| [out] | y | the output vector |
|
inline |
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blocks
| void Opm::gpuistl::GpuSparseMatrixGeneric< T >::preprocessSpMV | ( | ) |
Preprocess SpMV operation to optimize for sparsity pattern.
This function preprocesses the sparsity pattern of the matrix to optimize for the SpMV operation.
|
virtual |
umv computes y=Ax+y
| [in] | x | the vector to multiply with A |
| [in,out] | y | the vector to add and store the output in |
| void Opm::gpuistl::GpuSparseMatrixGeneric< T >::updateNonzeroValues | ( | const GpuSparseMatrixGeneric< T > & | matrix | ) |
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
| matrix | the matrix to extract the non-zero values from |
| void Opm::gpuistl::GpuSparseMatrixGeneric< T >::updateNonzeroValues | ( | const MatrixType & | matrix, |
| bool | copyNonZeroElementsDirectly = false ) |
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
| matrix | the matrix to extract the non-zero values from |
| copyNonZeroElementsDirectly | if 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. |
| MatrixType | is assumed to be a Dune::BCRSMatrix compatible matrix. |
|
virtual |
umv computes y=alpha * Ax + y
| [in] | alpha | The scaling factor for the matrix-vector product |
| [in] | x | the vector to multiply with A |
| [in,out] | y | the vector to add and store the output in |
|
staticconstexpr |
Maximum block size supported by this implementation.
This constant defines an upper bound on the block size to ensure reasonable compilation times. While this class itself could support larger values, functions that call dispatchOnBlocksize() might have limitations. This value can be increased if needed, but will increase compilation time due to template instantiations.