20#ifndef OPM_ROCSPARSESOLVER_BACKEND_HEADER_INCLUDED
21#define OPM_ROCSPARSESOLVER_BACKEND_HEADER_INCLUDED
25#include <opm/simulators/linalg/gpubridge/GpuResult.hpp>
26#include <opm/simulators/linalg/gpubridge/GpuSolver.hpp>
27#include <opm/simulators/linalg/gpubridge/WellContributions.hpp>
29#include <opm/simulators/linalg/gpubridge/rocm/rocsparsePreconditioner.hpp>
31#include <rocblas/rocblas.h>
32#include <rocsparse/rocsparse.h>
34#include <hip/hip_version.h>
36namespace Opm::Accelerator {
39template<
class Scalar,
unsigned int block_size>
48 using Base::verbosity;
49 using Base::platformID;
52 using Base::tolerance;
53 using Base::initialized;
58 bool useJacMatrix =
false;
60 bool analysis_done =
false;
61 std::shared_ptr<BlockedMatrix<Scalar>> mat{};
62 std::shared_ptr<BlockedMatrix<Scalar>> jacMat{};
64 rocsparse_direction dir = rocsparse_direction_row;
65 rocsparse_operation operation = rocsparse_operation_none;
66 rocsparse_handle handle;
67 rocblas_handle blas_handle;
68 rocsparse_mat_descr descr_A;
69#if HIP_VERSION >= 50400000
70 rocsparse_mat_info spmv_info;
74 rocsparse_int *d_Arows, *d_Acols;
76 Scalar *d_x, *d_b, *d_r, *d_rw, *d_p;
77 Scalar *d_pw, *d_s, *d_t, *d_v;
81 std::unique_ptr<rocsparsePreconditioner<Scalar, block_size> > prec;
96 void copy_system_to_gpu(Scalar* b);
100 void update_system_on_gpu(Scalar* vals, Scalar* b);
104 bool analyze_matrix();
108 bool create_preconditioner();
126 unsigned int platformID,
127 unsigned int deviceID,
128 std::string linsolver);
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition BlockedMatrix.hpp:29
This class is based on InverseOperatorResult struct from dune/istl/solver.hh It is needed to prevent ...
Definition GpuResult.hpp:31
GpuSolver(int linear_solver_verbosity, int max_it, Scalar tolerance_)
Construct a GpuSolver.
Definition GpuSolver.hpp:73
~rocsparseSolverBackend()
Destroy a openclSolver, and free memory.
Definition rocsparseSolverBackend.cpp:113
rocsparseSolverBackend(int linear_solver_verbosity, int maxit, Scalar tolerance, unsigned int platformID, unsigned int deviceID, std::string linsolver)
Construct a rocsparseSolver.
Definition rocsparseSolverBackend.cpp:62
void get_result(Scalar *x) override
Get result after linear solve, and peform postprocessing if necessary.
Definition rocsparseSolverBackend.cpp:638
rocsparseSolverBackend(int linear_solver_verbosity, int maxit, Scalar tolerance, bool opencl_ilu_reorder)
For the CPR coarse solver.
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51