22#ifndef OPM_ISTLSOLVER_HEADER_INCLUDED
23#define OPM_ISTLSOLVER_HEADER_INCLUDED
25#include <dune/istl/owneroverlapcopy.hh>
26#include <dune/istl/solver.hh>
28#include <opm/common/CriticalError.hpp>
29#include <opm/common/ErrorMacros.hpp>
30#include <opm/common/Exceptions.hpp>
31#include <opm/common/TimingMacros.hpp>
33#include <opm/grid/utility/ElementChunks.hpp>
39#include <opm/simulators/flow/BlackoilModelParameters.hpp>
42#include <opm/simulators/linalg/ExtractParallelGridInformationToISTL.hpp>
43#include <opm/simulators/linalg/FlowLinearSolverParameters.hpp>
44#include <opm/simulators/linalg/matrixblock.hh>
46#include <opm/simulators/linalg/PreconditionerWithUpdate.hpp>
47#include <opm/simulators/linalg/WellOperators.hpp>
48#include <opm/simulators/linalg/WriteSystemMatrixHelper.hpp>
49#include <opm/simulators/linalg/findOverlapRowsAndColumns.hpp>
50#include <opm/simulators/linalg/getQuasiImpesWeights.hpp>
51#include <opm/simulators/linalg/setupPropertyTree.hpp>
52#include <opm/simulators/linalg/AbstractISTLSolver.hpp>
53#include <opm/simulators/linalg/printlinearsolverparameter.hpp>
65namespace Opm::Properties {
69 using InheritsFrom = std::tuple<FlowIstlSolverParams>;
73template <
class TypeTag,
class MyTypeTag>
78template<
class TypeTag>
99template<
class Matrix,
class Vector,
class Comm>
100struct FlexibleSolverInfo
102 using AbstractSolverType = Dune::InverseOperator<Vector, Vector>;
103 using AbstractOperatorType = Dune::AssembledLinearOperator<Matrix, Vector, Vector>;
106 void create(
const Matrix& matrix,
109 std::size_t pressureIndex,
114 std::unique_ptr<AbstractSolverType> solver_;
115 std::unique_ptr<AbstractOperatorType> op_;
116 std::unique_ptr<LinearOperatorExtra<Vector,Vector>> wellOperator_;
117 AbstractPreconditionerType* pre_ =
nullptr;
118 std::size_t interiorCellNum_ = 0;
124void copyParValues(std::any& parallelInformation, std::size_t size,
125 Dune::OwnerOverlapCopyCommunication<int,int>& comm);
130template<
class Matrix>
131void makeOverlapRowsInvalid(Matrix& matrix,
136template<
class Matrix,
class Gr
id>
137std::unique_ptr<Matrix> blockJacobiAdjacency(
const Grid& grid,
139 std::size_t nonzeroes,
147 template <
class TypeTag>
158 using Matrix =
typename SparseMatrixAdapter::IstlMatrix;
161 using AbstractSolverType = Dune::InverseOperator<Vector, Vector>;
162 using AbstractOperatorType = Dune::AssembledLinearOperator<Matrix, Vector, Vector>;
171 constexpr static bool isIncompatibleWithCprw = enablePolymerMolarWeight;
174 using CommunicationType = Dune::OwnerOverlapCopyCommunication<int,int>;
176 using CommunicationType = Dune::Communication<int>;
180 using AssembledLinearOperatorType = Dune::AssembledLinearOperator< Matrix, Vector, Vector >;
182 static void registerParameters()
184 FlowLinearSolverParameters::registerParameters();
197 : simulator_(simulator),
209 : simulator_(simulator),
214 parameters_.resize(1);
215 parameters_[0].init(simulator_.vanguard().eclState().getSimulationConfig().useCPR());
223 if (isIncompatibleWithCprw) {
225 if (parameters_[0].linsolver_ ==
"cprw" || parameters_[0].linsolver_ ==
"hybrid") {
227 "The polymer injectivity model is incompatible with the CPRW linear solver.\n"
228 "Choose a different option, for example --linear-solver=ilu0");
232 if (parameters_[0].linsolver_ ==
"hybrid") {
240 FlowLinearSolverParameters
para;
242 para.linsolver_ =
"cprw";
243 parameters_.push_back(
para);
245 Parameters::IsSet<Parameters::LinearSolverMaxIter>(),
246 Parameters::IsSet<Parameters::LinearSolverReduction>()));
249 FlowLinearSolverParameters
para;
251 para.linsolver_ =
"ilu0";
252 parameters_.push_back(
para);
254 Parameters::IsSet<Parameters::LinearSolverMaxIter>(),
255 Parameters::IsSet<Parameters::LinearSolverReduction>()));
259 assert(parameters_.size() == 1);
263 if (parameters_[0].is_nldd_local_solver_) {
265 Parameters::IsSet<Parameters::NlddLocalLinearSolverMaxIter>(),
266 Parameters::IsSet<Parameters::NlddLocalLinearSolverReduction>()));
270 Parameters::IsSet<Parameters::LinearSolverMaxIter>(),
271 Parameters::IsSet<Parameters::LinearSolverReduction>()));
274 flexibleSolver_.resize(prm_.size());
276 const bool on_io_rank = (simulator_.gridView().comm().rank() == 0);
278 comm_.reset(
new CommunicationType( simulator_.vanguard().grid().comm() ) );
285 ElementMapper
elemMapper(simulator_.vanguard().gridView(), Dune::mcmgElementLayout());
286 detail::findOverlapAndInterior(simulator_.vanguard().grid(),
elemMapper, overlapRows_, interiorRows_);
287 useWellConn_ = Parameters::Get<Parameters::MatrixAddWellContributions>();
288 const bool ownersFirst = Parameters::Get<Parameters::OwnerCellsFirst>();
290 const std::string
msg =
"The linear solver no longer supports --owner-cells-first=false.";
297 const int interiorCellNum_ = detail::numMatrixRowsToUseInSolver(simulator_.vanguard().grid(),
true);
298 for (
auto&
f : flexibleSolver_) {
299 f.interiorCellNum_ = interiorCellNum_;
304 const std::size_t size = simulator_.vanguard().grid().leafGridView().size(0);
305 detail::copyParValues(parallelInformation_, size, *comm_);
310 detail::printLinearSolverParameters(parameters_, activeSolverNum_, prm_, simulator_.gridView().comm());
312 element_chunks_ = std::make_unique<ElementChunksType>(simulator_.vanguard().gridView(), Dune::Partitions::all,
ThreadManager::maxThreads());
322 if (
num >
static_cast<int>(prm_.size()) - 1) {
323 OPM_THROW(std::logic_error,
"Solver number " + std::to_string(
num) +
" not available.");
325 activeSolverNum_ =
num;
326 if (simulator_.gridView().comm().rank() == 0) {
327 OpmLog::debug(
"Active solver = " + std::to_string(activeSolverNum_)
328 +
" (" + parameters_[activeSolverNum_].linsolver_ +
")");
334 return flexibleSolver_.size();
337 void initPrepare(
const Matrix& M, Vector&
b)
339 const bool firstcall = (matrix_ ==
nullptr);
346 matrix_ =
const_cast<Matrix*
>(&M);
348 useWellConn_ = Parameters::Get<Parameters::MatrixAddWellContributions>();
352 if ( &M != matrix_ ) {
354 "Matrix objects are expected to be reused when reassembling!");
361 std::string type = prm_[activeSolverNum_].template
get<std::string>(
"preconditioner.type",
"paroverilu0");
362 std::transform(type.begin(), type.end(), type.begin(),
::tolower);
363 if (isParallel() && type !=
"paroverilu0") {
364 detail::makeOverlapRowsInvalid(getMatrix(), overlapRows_);
368 void prepare(
const SparseMatrixAdapter& M, Vector&
b)
override
379 prepareFlexibleSolver();
403 void resetSolveCount() {
412 const int verbosity = prm_[activeSolverNum_].get(
"verbosity", 0);
415 Helper::writeSystem(simulator_,
422 Dune::InverseOperatorResult
result;
425 assert(flexibleSolver_[activeSolverNum_].solver_);
426 flexibleSolver_[activeSolverNum_].solver_->apply(x, *rhs_,
result);
429 iterations_ =
result.iterations;
432 return checkConvergence(
result);
448 const CommunicationType*
comm()
const override {
return comm_.get(); }
450 void setDomainIndex(
const int index)
452 domainIndex_ = index;
455 bool isNlddLocalSolver()
const
457 return parameters_[activeSolverNum_].is_nldd_local_solver_;
462 using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
465 bool checkConvergence(
const Dune::InverseOperatorResult&
result)
const
470 bool isParallel()
const {
472 return !forceSerial_ && comm_->communicator().size() > 1;
478 void prepareFlexibleSolver()
483 if (isNlddLocalSolver()) {
484 auto wellOp = std::make_unique<DomainWellModelAsLinearOperator<WellModel, Vector, Vector>>(simulator_.problem().wellModel());
485 wellOp->setDomainIndex(domainIndex_);
486 flexibleSolver_[activeSolverNum_].wellOperator_ = std::move(
wellOp);
489 auto wellOp = std::make_unique<WellModelOperator>(simulator_.problem().wellModel());
490 flexibleSolver_[activeSolverNum_].wellOperator_ = std::move(
wellOp);
493 std::function<Vector()>
weightCalculator = this->getWeightsCalculator(prm_[activeSolverNum_], getMatrix(), pressureIndex);
495 flexibleSolver_[activeSolverNum_].create(getMatrix(),
497 prm_[activeSolverNum_],
506 flexibleSolver_[activeSolverNum_].pre_->update();
517 if (flexibleSolver_.empty()) {
520 if (!flexibleSolver_[activeSolverNum_].solver_) {
524 if (flexibleSolver_[activeSolverNum_].pre_->hasPerfectUpdate()) {
530 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 0) {
534 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 1) {
536 const int newton_iteration = this->simulator_.model().newtonMethod().numIterations();
539 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 2) {
543 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 3) {
547 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 4) {
549 const int step = this->parameters_[activeSolverNum_].cpr_reuse_interval_;
550 const bool create = ((solveCount_ % step) == 0);
554 const bool on_io_rank = (simulator_.gridView().comm().rank() == 0);
555 std::string
msg =
"Invalid value: " + std::to_string(this->parameters_[activeSolverNum_].cpr_reuse_setup_)
556 +
" for --cpr-reuse-setup parameter, run with --help to see allowed values.";
560 throw std::runtime_error(
msg);
569 std::function<Vector()> getWeightsCalculator(
const PropertyTree& prm,
570 const Matrix& matrix,
575 using namespace std::string_literals;
584 const auto weightsType = prm.
get(
"preconditioner.weight_type"s,
"quasiimpes"s);
590 return Amg::getQuasiImpesWeights<Matrix, Vector>(matrix,
599 Vector weights(rhs_->size());
600 ElementContext elemCtx(simulator_);
614 Vector weights(rhs_->size());
615 ElementContext elemCtx(simulator_);
628 "not implemented for cpr."
629 " Please use quasiimpes, trueimpes or trueimpesanalytic.");
641 const Matrix& getMatrix()
const
646 const Simulator& simulator_;
647 mutable int iterations_;
648 mutable int solveCount_;
649 std::any parallelInformation_;
655 int activeSolverNum_ = 0;
656 std::vector<detail::FlexibleSolverInfo<Matrix,Vector,CommunicationType>> flexibleSolver_;
657 std::vector<int> overlapRows_;
658 std::vector<int> interiorRows_;
660 int domainIndex_ = -1;
664 std::vector<FlowLinearSolverParameters> parameters_;
665 bool forceSerial_ =
false;
666 std::vector<PropertyTree> prm_;
668 std::shared_ptr< CommunicationType > comm_;
669 std::unique_ptr<ElementChunksType> element_chunks_;
Helper class for grid instantiation of ECL file-format using problems.
Interface class adding the update() method to the preconditioner interface.
Definition PreconditionerWithUpdate.hpp:32
Abstract interface for ISTL solvers.
Definition AbstractISTLSolver.hpp:45
static bool checkConvergence(const Dune::InverseOperatorResult &result, const FlowLinearSolverParameters ¶meters)
Check the convergence of the linear solver.
Definition AbstractISTLSolver.hpp:194
This class solves the fully implicit black-oil system by solving the reduced system (after eliminatin...
Definition ISTLSolver.hpp:149
void getResidual(Vector &b) const override
Get the residual vector.
Definition ISTLSolver.hpp:389
ISTLSolver(const Simulator &simulator, const FlowLinearSolverParameters ¶meters, bool forceSerial=false)
Construct a system solver.
Definition ISTLSolver.hpp:194
int iterations() const override
Solve the system of linear equations Ax = b, with A being the combined derivative matrix of the resid...
Definition ISTLSolver.hpp:443
void setActiveSolver(const int num) override
Set the active solver by its index.
Definition ISTLSolver.hpp:320
void setMatrix(const SparseMatrixAdapter &) override
Set the matrix for the solver.
Definition ISTLSolver.hpp:394
void prepare(const Matrix &M, Vector &b) override
Prepare the solver with the given matrix and right-hand side vector.
Definition ISTLSolver.hpp:373
int numAvailableSolvers() const override
Get the number of available solvers.
Definition ISTLSolver.hpp:332
void prepare(const SparseMatrixAdapter &M, Vector &b) override
Prepare the solver with the given sparse matrix and right-hand side vector.
Definition ISTLSolver.hpp:368
int getSolveCount() const override
Get the count of how many times the solver has been called.
Definition ISTLSolver.hpp:399
bool shouldCreateSolver() const
Return true if we should (re)create the whole solver, instead of just calling update() on the precond...
Definition ISTLSolver.hpp:513
const std::any & parallelInformation() const
Definition ISTLSolver.hpp:446
void eraseMatrix() override
Signals that the memory for the matrix internally in the solver could be erased.
Definition ISTLSolver.hpp:316
ISTLSolver(const Simulator &simulator)
Construct a system solver.
Definition ISTLSolver.hpp:208
bool solve(Vector &x) override
Solve the system of equations Ax = b.
Definition ISTLSolver.hpp:407
void setResidual(Vector &) override
Set the residual vector.
Definition ISTLSolver.hpp:384
const CommunicationType * comm() const override
Get the communication object used by the solver.
Definition ISTLSolver.hpp:448
A sparse matrix interface backend for BCRSMatrix from dune-istl.
Definition istlsparsematrixadapter.hh:43
Definition matrixblock.hh:229
Hierarchical collection of key/value pairs.
Definition PropertyTree.hpp:39
T get(const std::string &key) const
Retrieve property value given hierarchical property key.
Definition PropertyTree.cpp:59
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition threadmanager.hpp:66
Definition WellOperators.hpp:70
Declare the properties used by the infrastructure code of the finite volume discretizations.
A sparse matrix interface backend for BCRSMatrix from dune-istl.
Defines the common properties required by the porous medium multi-phase models.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233
PropertyTree setupPropertyTree(FlowLinearSolverParameters p, bool linearSolverMaxIterSet, bool linearSolverReductionSet)
Set up a property tree intended for FlexibleSolver by either reading the tree from a JSON file or cre...
Definition setupPropertyTree.cpp:182
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
This class carries all parameters for the NewtonIterationBlackoilInterleaved class.
Definition FlowLinearSolverParameters.hpp:98
The class that allows to manipulate sparse matrices.
Definition linalgproperties.hh:50
Definition ISTLSolver.hpp:68
Definition FlowBaseProblemProperties.hpp:97