20#ifndef OPM_AMGCLSOLVER_BACKEND_HEADER_INCLUDED
21#define OPM_AMGCLSOLVER_BACKEND_HEADER_INCLUDED
23#include <opm/simulators/linalg/gpubridge/GpuResult.hpp>
24#include <opm/simulators/linalg/gpubridge/GpuSolver.hpp>
25#include <opm/simulators/linalg/gpubridge/WellContributions.hpp>
27#include <boost/property_tree/ptree.hpp>
29#include <amgcl/amg.hpp>
30#include <amgcl/backend/builtin.hpp>
31#include <amgcl/adapter/crs_tuple.hpp>
32#include <amgcl/make_block_solver.hpp>
33#include <amgcl/relaxation/as_preconditioner.hpp>
34#include <amgcl/relaxation/ilu0.hpp>
35#include <amgcl/solver/bicgstab.hpp>
36#include <amgcl/preconditioner/runtime.hpp>
37#include <amgcl/value_type/static_matrix.hpp>
44namespace Opm::Accelerator {
48template<
class Scalar,
unsigned int block_size>
57 using Base::verbosity;
58 using Base::platformID;
61 using Base::tolerance;
62 using Base::initialized;
64 using dmat_type = amgcl::static_matrix<Scalar, block_size, block_size>;
65 using dvec_type = amgcl::static_matrix<Scalar, block_size, 1>;
66 using CPU_Backend = std::conditional_t<block_size == 1,
67 amgcl::backend::builtin<Scalar>,
68 amgcl::backend::builtin<dmat_type>>;
70 using CPU_Solver = amgcl::make_solver<amgcl::runtime::preconditioner<CPU_Backend>,
71 amgcl::runtime::solver::wrapper<CPU_Backend>>;
75 enum Amgcl_backend_type {
82 std::vector<unsigned> A_rows, A_cols;
83 std::vector<Scalar> A_vals, rhs;
84 std::vector<Scalar> x;
85 std::once_flag print_info;
86 Amgcl_backend_type backend_type = cpu;
103 void initialize(
int Nb,
int nnzbs);
108 void convert_sparsity_pattern(
const int *rows,
const int *cols);
113 void convert_data(
const Scalar*
vals,
const int* rows);
128 Scalar tolerance,
unsigned int platformID,
129 unsigned int deviceID);
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
This class serves to simplify choosing between different backend solvers, such as cusparseSolver and ...
Definition GpuSolver.hpp:46
This class does not implement a solver, but converts the BCSR format to normal CSR and uses amgcl for...
Definition amgclSolverBackend.hpp:50
void get_result(Scalar *x) override
Get result after linear solve, and peform postprocessing if necessary.
Definition amgclSolverBackend.cpp:401
~amgclSolverBackend() override
Destroy a openclSolver, and free memory.
Definition amgclSolverBackend.cpp:64
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51
Definition PropertyTree.hpp:31
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242