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>
31 class WellContributions;
33 namespace Accelerator {
34 enum class SolverStatus {
36 BDA_SOLVER_ANALYSIS_FAILED,
37 BDA_SOLVER_CREATE_PRECONDITIONER_FAILED,
38 BDA_SOLVER_UNKNOWN_ERROR
43 template <
unsigned int block_size>
58 double tolerance = 1e-2;
60 std::string bitstream;
67 unsigned int platformID = 0;
68 unsigned int deviceID = 0;
70 bool initialized =
false;
80 BdaSolver(
int linear_solver_verbosity,
int max_it,
double tolerance_) : verbosity(linear_solver_verbosity), maxit(max_it), tolerance(tolerance_) {};
81 BdaSolver(
int linear_solver_verbosity,
int max_it,
double tolerance_,
unsigned int deviceID_) : verbosity(linear_solver_verbosity), maxit(max_it), tolerance(tolerance_), deviceID(deviceID_) {};
82 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_) {};
83 BdaSolver(std::string fpga_bitstream,
int linear_solver_verbosity,
int max_it,
double tolerance_) : verbosity(linear_solver_verbosity), maxit(max_it), tolerance(tolerance_), bitstream(fpga_bitstream) {};
89 virtual SolverStatus
solve_system(std::shared_ptr<BlockedMatrix> matrix,
double *b,
92 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:45
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:86
BdaSolver(int linear_solver_verbosity, int max_it, double tolerance_)
Construct a BdaSolver, can be cusparseSolver, openclSolver, fpgaSolver.
Definition: BdaSolver.hpp:80
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition: WellContributions.hpp:53
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27