opm-simulators
Loading...
Searching...
No Matches
Opm::Accelerator::rocsparseBILU0< Scalar, block_size > Class Template Reference

This class implements a Blocked ILU0 preconditioner The decomposition is done on GPU, using exact decomposition, or ChowPatel decomposition The preconditioner is applied via two exact triangular solves. More...

#include <rocsparseBILU0.hpp>

Inheritance diagram for Opm::Accelerator::rocsparseBILU0< Scalar, block_size >:
Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size > Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >

Public Member Functions

 rocsparseBILU0 (int verbosity_)
bool initialize (std::shared_ptr< BlockedMatrix< Scalar > > matrix, std::shared_ptr< BlockedMatrix< Scalar > > jacMatrix, rocsparse_int *d_Arows, rocsparse_int *d_Acols) override
 Initialize GPU and allocate memory.
bool analyze_matrix ()
 Analysis, extract parallelism if specified.
bool analyze_matrix (BlockedMatrix< Scalar > *mat) override
 Analysis, extract parallelism if specified.
bool analyze_matrix (BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
 Analysis, extract parallelism if specified.
bool create_preconditioner (BlockedMatrix< Scalar > *mat) override
 ILU decomposition.
bool create_preconditioner (BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
 ILU decomposition.
void apply (const Scalar &y, Scalar &x, WellContributions< Scalar > &wellContribs) override
 Apply preconditioner, x = prec(y) via Lz = y and Ux = z.
void copy_system_to_gpu (Scalar *mVals) override
 Copy matrix A values to GPU.
void copy_values_to_gpu (Scalar *mVals, int *mRows, int *mCols, bool reuse)
 Copy matrix A values to GPU.
void update_system_on_gpu (Scalar *, Scalar *b) override
 Update GPU values after a new assembly is done.
Public Member Functions inherited from Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >
void set_matrix_analysis (rocsparse_mat_descr descr_L, rocsparse_mat_descr descr_U)
void set_context (rocsparse_handle handle, rocblas_handle blas_handle, rocsparse_direction dir, rocsparse_operation operation, hipStream_t stream)
void setJacMat (const BlockedMatrix< Scalar > &jacMat)
Public Member Functions inherited from Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >
virtual void apply (const ApplyScalar &y, ApplyScalar &x, WellContributions< Scalar > &wellContribs)=0

Additional Inherited Members

Static Public Member Functions inherited from Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >
static std::unique_ptr< rocsparsePreconditioner< Scalar, block_size > > create (PreconditionerType type, int verbosity)
Static Public Member Functions inherited from Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >
static std::unique_ptr< Preconditioner > create (PreconditionerType type, bool opencl_ilu_parallel, int verbosity)
Public Attributes inherited from Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >
int nnzbs_prec = 0
bool useJacMatrix = false
std::shared_ptr< BlockedMatrix< Scalar > > jacMat {}
Protected Member Functions inherited from Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >
 rocsparsePreconditioner (int verbosity_)
Protected Member Functions inherited from Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >
 Preconditioner (int verbosity_)
Protected Attributes inherited from Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >
rocsparse_handle handle
rocblas_handle blas_handle
rocsparse_direction dir = rocsparse_direction_row
rocsparse_operation operation = rocsparse_operation_none
rocsparse_mat_descr descr_L
rocsparse_mat_descr descr_U
hipStream_t stream
Protected Attributes inherited from Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >
int N = 0
int Nb = 0
int nnz = 0
int nnzb = 0
int verbosity = 0

Detailed Description

template<class Scalar, unsigned int block_size>
class Opm::Accelerator::rocsparseBILU0< Scalar, block_size >

This class implements a Blocked ILU0 preconditioner The decomposition is done on GPU, using exact decomposition, or ChowPatel decomposition The preconditioner is applied via two exact triangular solves.

Member Function Documentation

◆ analyze_matrix() [1/2]

template<class Scalar, unsigned int block_size>
bool Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::analyze_matrix ( BlockedMatrix< Scalar > * mat)
overridevirtual

Analysis, extract parallelism if specified.

Parameters
[in]matmatrix A

Implements Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >.

◆ analyze_matrix() [2/2]

template<class Scalar, unsigned int block_size>
bool Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::analyze_matrix ( BlockedMatrix< Scalar > * mat,
BlockedMatrix< Scalar > * jacMat )
overridevirtual

Analysis, extract parallelism if specified.

Parameters
[in]matmatrix A
[in]jacMatmatrix for preconditioner, analyze this as well

Implements Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >.

◆ apply()

template<class Scalar, unsigned int block_size>
void Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::apply ( const Scalar & y,
Scalar & x,
WellContributions< Scalar > & wellContribs )
override

Apply preconditioner, x = prec(y) via Lz = y and Ux = z.

Parameters
[in]yInput y vector
[out]xOutput x vector
wellContribsWell contributions

◆ copy_system_to_gpu()

template<class Scalar, unsigned int block_size>
void Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::copy_system_to_gpu ( Scalar * mVals)
overridevirtual

Copy matrix A values to GPU.

Parameters
[in]mValsInput values

Implements Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >.

◆ copy_values_to_gpu()

template<class Scalar, unsigned int block_size>
void Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::copy_values_to_gpu ( Scalar * mVals,
int * mRows,
int * mCols,
bool reuse )

Copy matrix A values to GPU.

Parameters
[in]mValsInput values
[in]mRowsArray of matrix row indices
[in]mColsArray of matrix column indices
[in]reuseTrue to reuse old matrix

◆ create_preconditioner() [1/2]

template<class Scalar, unsigned int block_size>
bool Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::create_preconditioner ( BlockedMatrix< Scalar > * mat)
overridevirtual

ILU decomposition.

Parameters
[in]matmatrix A to decompose

Implements Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >.

◆ create_preconditioner() [2/2]

template<class Scalar, unsigned int block_size>
bool Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::create_preconditioner ( BlockedMatrix< Scalar > * mat,
BlockedMatrix< Scalar > * jacMat )
overridevirtual

ILU decomposition.

Parameters
[in]matmatrix A
[in]jacMatmatrix for preconditioner, decompose this one if used

Implements Opm::Accelerator::Preconditioner< Scalar, block_size, ApplyScalar >.

◆ initialize()

template<class Scalar, unsigned int block_size>
bool Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::initialize ( std::shared_ptr< BlockedMatrix< Scalar > > matrix,
std::shared_ptr< BlockedMatrix< Scalar > > jacMatrix,
rocsparse_int * d_Arows,
rocsparse_int * d_Acols )
overridevirtual

Initialize GPU and allocate memory.

Parameters
[in]matrixmatrix A
[in]jacMatrixmatrix for preconditioner
[in]d_ArowsArray of row indices
[in]d_AcolsArray of column indices

Implements Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >.

◆ update_system_on_gpu()

template<class Scalar, unsigned int block_size>
void Opm::Accelerator::rocsparseBILU0< Scalar, block_size >::update_system_on_gpu ( Scalar * ,
Scalar * b )
overridevirtual

Update GPU values after a new assembly is done.

Parameters
[in]bNew b vector

Implements Opm::Accelerator::rocsparsePreconditioner< Scalar, block_size >.


The documentation for this class was generated from the following files:
  • opm/simulators/linalg/gpubridge/rocm/rocsparseBILU0.hpp
  • opm/simulators/linalg/gpubridge/rocm/rocsparseBILU0.cpp