20#ifndef OPM_BDASOLVER_BACKEND_HEADER_INCLUDED
21#define OPM_BDASOLVER_BACKEND_HEADER_INCLUDED
24#include <opm/simulators/linalg/bda/BdaResult.hpp>
25#include <opm/simulators/linalg/bda/BlockedMatrix.hpp>
32class WellContributions;
34namespace Accelerator {
35 enum class SolverStatus {
37 BDA_SOLVER_ANALYSIS_FAILED,
38 BDA_SOLVER_CREATE_PRECONDITIONER_FAILED,
39 BDA_SOLVER_UNKNOWN_ERROR
44 template <
unsigned int block_size>
59 double tolerance = 1e-2;
66 unsigned int platformID = 0;
67 unsigned int deviceID = 0;
69 bool initialized =
false;
78 BdaSolver(
int linear_solver_verbosity,
int max_it,
double tolerance_) : verbosity(linear_solver_verbosity), maxit(max_it), tolerance(tolerance_) {};
79 BdaSolver(
int linear_solver_verbosity,
int max_it,
double tolerance_,
unsigned int deviceID_) : verbosity(linear_solver_verbosity), maxit(max_it), tolerance(tolerance_), deviceID(deviceID_) {};
80 BdaSolver(
int linear_solver_verbosity,
int max_it,
double tolerance_,
unsigned int platformID_,
unsigned int deviceID_) : verbosity(linear_solver_verbosity), maxit(max_it), tolerance(tolerance_), platformID(platformID_), deviceID(deviceID_) {};
86 virtual SolverStatus
solve_system(std::shared_ptr<BlockedMatrix> matrix,
double *b,
89 virtual void get_result(
double *x) = 0;
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
virtual SolverStatus solve_system(std::shared_ptr< BlockedMatrix > matrix, double *b, std::shared_ptr< BlockedMatrix > jacMatrix, WellContributions &wellContribs, BdaResult &res)=0
Define as pure virtual functions, so derivedclass must implement them.
virtual ~BdaSolver()
Define virtual destructor, so that the derivedclass destructor will be called.
Definition BdaSolver.hpp:83
BdaSolver(int linear_solver_verbosity, int max_it, double tolerance_)
Construct a BdaSolver.
Definition BdaSolver.hpp:78
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:52
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27