21#ifndef OPM_NON_LINEAR_SOLVER_EBOS_HPP
22#define OPM_NON_LINEAR_SOLVER_EBOS_HPP
24#include <opm/simulators/timestepping/SimulatorReport.hpp>
25#include <opm/common/ErrorMacros.hpp>
26#include <opm/simulators/timestepping/SimulatorTimerInterface.hpp>
28#include <opm/models/utils/parametersystem.hh>
29#include <opm/models/utils/propertysystem.hh>
30#include <opm/models/utils/basicproperties.hh>
31#include <opm/models/nonlinear/newtonmethodproperties.hh>
32#include <opm/common/Exceptions.hpp>
34#include <dune/common/fmatrix.hh>
35#include <dune/istl/bcrsmatrix.hh>
38namespace Opm::Properties {
44template<
class TypeTag,
class MyTypeTag>
46 using type = UndefinedProperty;
52template<
class TypeTag,
class MyTypeTag>
54 using type = UndefinedProperty;
56template<
class TypeTag,
class MyTypeTag>
58 using type = UndefinedProperty;
61template<
class TypeTag>
63 using type = GetPropType<TypeTag, Scalar>;
64 static constexpr type value = 0.5;
66template<
class TypeTag>
67struct NewtonMaxIterations<TypeTag, TTag::FlowNonLinearSolver> {
68 static constexpr int value = 20;
70template<
class TypeTag>
72 static constexpr int value = 1;
74template<
class TypeTag>
76 static constexpr auto value =
"dampen";
86enum class NonlinearRelaxType {
94void detectOscillations(
const std::vector<std::vector<double>>& residualHistory,
95 const int it,
const int numPhases,
const double relaxRelTol,
96 bool& oscillate,
bool& stagnate);
100template <
class BVector>
101void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
102 const double omega, NonlinearRelaxType relaxType);
108 template <
class TypeTag,
class PhysicalModel>
111 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
117 NonlinearRelaxType relaxType_;
119 double relaxIncrement_;
130 relaxMax_ = EWOMS_GET_PARAM(TypeTag, Scalar, NewtonMaxRelax);
131 maxIter_ = EWOMS_GET_PARAM(TypeTag,
int, NewtonMaxIterations);
132 minIter_ = EWOMS_GET_PARAM(TypeTag,
int, NewtonMinIterations);
134 const auto& relaxationTypeString = EWOMS_GET_PARAM(TypeTag, std::string, NewtonRelaxationType);
135 if (relaxationTypeString ==
"dampen") {
136 relaxType_ = NonlinearRelaxType::Dampen;
137 }
else if (relaxationTypeString ==
"sor") {
138 relaxType_ = NonlinearRelaxType::SOR;
140 OPM_THROW(std::runtime_error,
141 "Unknown Relaxtion Type " + relaxationTypeString);
145 static void registerParameters()
147 EWOMS_REGISTER_PARAM(TypeTag, Scalar, NewtonMaxRelax,
"The maximum relaxation factor of a Newton iteration");
148 EWOMS_REGISTER_PARAM(TypeTag,
int, NewtonMaxIterations,
"The maximum number of Newton iterations per time step");
149 EWOMS_REGISTER_PARAM(TypeTag,
int, NewtonMinIterations,
"The minimum number of Newton iterations per time step");
150 EWOMS_REGISTER_PARAM(TypeTag, std::string, NewtonRelaxationType,
"The type of relaxation used by Newton method");
156 relaxType_ = NonlinearRelaxType::Dampen;
158 relaxIncrement_ = 0.1;
179 std::unique_ptr<PhysicalModel>
model)
181 , model_(std::move(
model))
183 , nonlinearIterations_(0)
184 , linearIterations_(0)
186 , nonlinearIterationsLast_(0)
187 , linearIterationsLast_(0)
188 , wellIterationsLast_(0)
191 OPM_THROW(std::logic_error,
"Must provide a non-null model argument for NonlinearSolver.");
203 report += model_->prepareStep(timer);
210 bool converged =
false;
218 auto iterReport = model_->nonlinearIteration(iteration, timer, *
this);
220 report += iterReport;
221 report.converged = iterReport.converged;
223 converged = report.converged;
229 failureReport_ = report;
230 failureReport_ += model_->failureReport();
234 while ( (!converged && (iteration <=
maxIter())) || (iteration <=
minIter()));
237 failureReport_ = report;
239 std::string msg =
"Solver convergence failure - Failed to complete a time step within " + std::to_string(
maxIter()) +
" iterations.";
240 OPM_THROW_NOLOG(TooManyIterations, msg);
244 report += model_->afterStep(timer);
245 report.converged =
true;
251 {
return failureReport_; }
255 {
return linearizations_; }
259 {
return nonlinearIterations_; }
263 {
return linearIterations_; }
267 {
return wellIterations_; }
271 {
return nonlinearIterationsLast_; }
275 {
return linearIterationsLast_; }
279 {
return wellIterationsLast_; }
281 std::vector<std::vector<double> >
282 computeFluidInPlace(
const std::vector<int>& fipnum)
const
283 {
return model_->computeFluidInPlace(fipnum); }
295 const int it,
bool& oscillate,
bool& stagnate)
const
297 detail::detectOscillations(residualHistory, it, model_->numPhases(),
298 this->relaxRelTol(), oscillate, stagnate);
304 template <
class BVector>
307 detail::stabilizeNonlinearUpdate(dx, dxOld, omega, this->
relaxType());
312 {
return param_.relaxMax_; }
316 {
return param_.relaxIncrement_; }
320 {
return param_.relaxType_; }
324 {
return param_.relaxRelTol_; }
328 {
return param_.maxIter_; }
332 {
return param_.minIter_; }
341 SolverParameters param_;
342 std::unique_ptr<PhysicalModel> model_;
344 int nonlinearIterations_;
345 int linearIterations_;
347 int nonlinearIterationsLast_;
348 int linearIterationsLast_;
349 int wellIterationsLast_;
A nonlinear solver class suitable for general fully-implicit models, as well as pressure,...
Definition: NonlinearSolverEbos.hpp:110
double relaxMax() const
The greatest relaxation factor (i.e. smallest factor) allowed.
Definition: NonlinearSolverEbos.hpp:311
int maxIter() const
The maximum number of nonlinear iterations allowed.
Definition: NonlinearSolverEbos.hpp:327
const SimulatorReportSingle & failureReport() const
return the statistics if the step() method failed
Definition: NonlinearSolverEbos.hpp:250
PhysicalModel & model()
Mutable reference to physical model.
Definition: NonlinearSolverEbos.hpp:290
const PhysicalModel & model() const
Reference to physical model.
Definition: NonlinearSolverEbos.hpp:286
int wellIterationsLastStep() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:278
void detectOscillations(const std::vector< std::vector< double > > &residualHistory, const int it, bool &oscillate, bool &stagnate) const
Detect oscillation or stagnation in a given residual history.
Definition: NonlinearSolverEbos.hpp:294
double relaxIncrement() const
The step-change size for the relaxation factor.
Definition: NonlinearSolverEbos.hpp:315
int linearIterationsLastStep() const
Number of linear solver iterations used in the last call to step().
Definition: NonlinearSolverEbos.hpp:274
int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:262
int linearizations() const
Number of linearizations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:254
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const double omega) const
Apply a stabilization to dx, depending on dxOld and relaxation parameters.
Definition: NonlinearSolverEbos.hpp:305
int minIter() const
The minimum number of nonlinear iterations allowed.
Definition: NonlinearSolverEbos.hpp:331
int nonlinearIterationsLastStep() const
Number of nonlinear solver iterations used in the last call to step().
Definition: NonlinearSolverEbos.hpp:270
double relaxRelTol() const
The relaxation relative tolerance.
Definition: NonlinearSolverEbos.hpp:323
int wellIterations() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:266
NonlinearSolverEbos(const SolverParameters ¶m, std::unique_ptr< PhysicalModel > model)
Construct solver for a given model.
Definition: NonlinearSolverEbos.hpp:178
void setParameters(const SolverParameters ¶m)
Set parameters to override those given at construction time.
Definition: NonlinearSolverEbos.hpp:335
NonlinearRelaxType relaxType() const
The relaxation type (Dampen or SOR).
Definition: NonlinearSolverEbos.hpp:319
int nonlinearIterations() const
Number of full nonlinear solver iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:258
Interface class for SimulatorTimer objects, to be improved.
Definition: SimulatorTimerInterface.hpp:34
virtual double currentStepLength() const =0
Current step length.
virtual double simulationTimeElapsed() const =0
Time elapsed since the start of the simulation until the beginning of the current time step [s].
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition: WellState.hpp:60
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27
Definition: NonlinearSolverEbos.hpp:116
Definition: NonlinearSolverEbos.hpp:45
Definition: NonlinearSolverEbos.hpp:53
Definition: NonlinearSolverEbos.hpp:57
Definition: NonlinearSolverEbos.hpp:41
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34