20#ifndef OPM_AMGCLSOLVER_BACKEND_HEADER_INCLUDED
21#define OPM_AMGCLSOLVER_BACKEND_HEADER_INCLUDED
23#include <opm/simulators/linalg/bda/BdaResult.hpp>
24#include <opm/simulators/linalg/bda/BdaSolver.hpp>
25#include <opm/simulators/linalg/bda/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>
50template <
unsigned int block_size>
59 using Base::verbosity;
60 using Base::platformID;
63 using Base::tolerance;
64 using Base::initialized;
66 typedef amgcl::static_matrix<double, block_size, block_size> dmat_type;
67 typedef amgcl::static_matrix<double, block_size, 1> dvec_type;
68 typedef amgcl::backend::builtin<dmat_type> CPU_Backend;
70 typedef amgcl::make_solver<amgcl::runtime::preconditioner<CPU_Backend>, amgcl::runtime::solver::wrapper<CPU_Backend> > CPU_Solver;
75 enum Amgcl_backend_type {
82 std::vector<unsigned> A_rows, A_cols;
83 std::vector<double> A_vals, rhs;
84 std::vector<double> 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(
int *rows,
int *cols);
113 void convert_data(
double *
vals,
int *rows);
139 SolverStatus solve_system(std::shared_ptr<BlockedMatrix> matrix,
double *
b,
This class is based on InverseOperatorResult struct from dune/istl/solver.hh It is needed to prevent ...
Definition BdaResult.hpp:31
This class serves to simplify choosing between different backend solvers, such as cusparseSolver and ...
Definition BdaSolver.hpp:46
This class does not implement a solver, but converts the BCSR format to normal CSR and uses amgcl for...
Definition amgclSolverBackend.hpp:52
~amgclSolverBackend()
Destroy a openclSolver, and free memory.
Definition amgclSolverBackend.cpp:67
void get_result(double *x) override
Get result after linear solve, and peform postprocessing if necessary.
Definition amgclSolverBackend.cpp:365
Definition AquiferInterface.hpp:35
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:52
Definition PropertyTree.hpp:29
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27