20#ifndef OPM_ROCSPARSEBILU0_HPP
21#define OPM_ROCSPARSEBILU0_HPP
23#include <opm/simulators/linalg/gpubridge/BlockedMatrix.hpp>
24#include <opm/simulators/linalg/gpubridge/WellContributions.hpp>
26#include <opm/simulators/linalg/gpubridge/rocm/rocsparsePreconditioner.hpp>
28#include <rocblas/rocblas.h>
29#include <rocsparse/rocsparse.h>
31#include <hip/hip_version.h>
33namespace Opm::Accelerator {
38template <
class Scalar,
unsigned int block_size>
47 using Base::verbosity;
53#if HIP_VERSION >= 50400000
58 Scalar *d_Mvals, *d_t;
61 std::size_t d_bufferSize_M=0, d_bufferSize_L=0, d_bufferSize_U=0, d_bufferSize=0;
107 void apply(
const Scalar&
y,
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition BlockedMatrix.hpp:29
This class implements a Blocked ILU0 preconditioner The decomposition is done on GPU,...
Definition rocsparseBILU0.hpp:40
bool create_preconditioner(BlockedMatrix< Scalar > *mat) override
ILU decomposition.
Definition rocsparseBILU0.cpp:328
void apply(const Scalar &y, Scalar &x, WellContributions< Scalar > &wellContribs) override
Apply preconditioner, x = prec(y) via Lz = y and Ux = z.
Definition rocsparseBILU0.cpp:428
void copy_values_to_gpu(Scalar *mVals, int *mRows, int *mCols, bool reuse)
Copy matrix A values to GPU.
Definition rocsparseBILU0.cpp:393
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.
Definition rocsparseBILU0.cpp:75
bool analyze_matrix()
Analysis, extract parallelism if specified.
Definition rocsparseBILU0.cpp:124
void copy_system_to_gpu(Scalar *mVals) override
Copy matrix A values to GPU.
Definition rocsparseBILU0.cpp:369
void update_system_on_gpu(Scalar *, Scalar *b) override
Update GPU values after a new assembly is done.
Definition rocsparseBILU0.cpp:406
Definition rocsparsePreconditioner.hpp:34
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240