21#ifndef OPM_NON_LINEAR_SOLVER_HPP
22#define OPM_NON_LINEAR_SOLVER_HPP
24#include <dune/common/fmatrix.hh>
25#include <dune/istl/bcrsmatrix.hh>
27#include <opm/common/ErrorMacros.hpp>
28#include <opm/common/Exceptions.hpp>
30#include <opm/models/nonlinear/newtonmethodparams.hpp>
31#include <opm/models/nonlinear/newtonmethodproperties.hh>
36#include <opm/simulators/timestepping/SimulatorReport.hpp>
37#include <opm/simulators/timestepping/SimulatorTimerInterface.hpp>
38#include <opm/simulators/timestepping/TimeStepControl.hpp>
42namespace Opm::Parameters {
54enum class NonlinearRelaxType {
63void detectOscillations(
const std::vector<std::vector<Scalar>>& residualHistory,
64 const int it,
const int numPhases,
const Scalar relaxRelTol,
65 const int minimumOscillatingPhases,
66 bool& oscillate,
bool& stagnate);
70template <
class BVector,
class Scalar>
71void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
72 const Scalar omega, NonlinearRelaxType relaxType);
78struct NonlinearSolverParameters
80 NonlinearRelaxType relaxType_;
82 Scalar relaxIncrement_;
85 NonlinearSolverParameters();
87 static void registerParameters();
94 template <
class TypeTag,
class PhysicalModel>
112 std::unique_ptr<PhysicalModel>
model)
114 , model_(std::move(
model))
116 , nonlinearIterations_(0)
117 , linearIterations_(0)
119 , nonlinearIterationsLast_(0)
120 , linearIterationsLast_(0)
121 , wellIterationsLast_(0)
124 OPM_THROW(std::logic_error,
"Must provide a non-null model argument for NonlinearSolver.");
136 report += model_->prepareStep(timer);
139 bool converged =
false;
147 auto iterReport = model_->nonlinearIteration(timer, *
this);
149 report += iterReport;
150 report.converged = iterReport.converged;
152 converged = report.converged;
157 failureReport_ = report;
158 failureReport_ += model_->failureReport();
162 while ( (!converged && (model_->simulator().problem().iterationContext().iteration() <= this->model().param().newton_max_iter_)) ||
163 (model_->simulator().problem().iterationContext().iteration() <= this->model().param().newton_min_iter_));
166 failureReport_ = report;
168 std::string msg =
"Solver convergence failure - Failed to complete a time step within ";
169 msg += std::to_string(
model().param().newton_max_iter_) +
" iterations.";
170 OPM_THROW_NOLOG(TooManyIterations, msg);
172 auto relativeChange = model_->relativeChange();
174 report.converged =
false;
175 report.time_step_rejected =
true;
176 failureReport_ = report;
178 std::string msg =
"Relative change in solution for time step was " + std::to_string(relativeChange);
179 msg +=
", which is larger than the tolerance accepted by the timestepping algorithm.";
180 OPM_THROW_NOLOG(TimeSteppingBreakdown, msg);
183 report.converged =
true;
189 {
return failureReport_; }
193 {
return linearizations_; }
197 {
return nonlinearIterations_; }
201 {
return linearIterations_; }
205 {
return wellIterations_; }
209 {
return nonlinearIterationsLast_; }
213 {
return linearIterationsLast_; }
217 {
return wellIterationsLast_; }
219 std::vector<std::vector<Scalar> >
220 computeFluidInPlace(
const std::vector<int>& fipnum)
const
221 {
return model_->computeFluidInPlace(fipnum); }
233 const int it,
bool& oscillate,
bool& stagnate)
const
235 detail::detectOscillations(residualHistory, it, model_->numPhases(),
236 this->relaxRelTol(), 1, oscillate, stagnate);
242 template <
class BVector>
245 detail::stabilizeNonlinearUpdate(dx, dxOld, omega, this->
relaxType());
250 {
return param_.relaxMax_; }
254 {
return param_.relaxIncrement_; }
258 {
return param_.relaxType_; }
262 {
return param_.relaxRelTol_; }
271 SolverParameters param_;
272 std::unique_ptr<PhysicalModel> model_;
274 int nonlinearIterations_;
275 int linearIterations_;
277 int nonlinearIterationsLast_;
278 int linearIterationsLast_;
279 int wellIterationsLast_;
Defines a type tags and some fundamental properties all models.
int linearIterationsLastStep() const
Number of linear solver iterations used in the last call to step().
Definition NonlinearSolver.hpp:212
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const Scalar omega) const
Apply a stabilization to dx, depending on dxOld and relaxation parameters.
Definition NonlinearSolver.hpp:243
Scalar relaxRelTol() const
The relaxation relative tolerance.
Definition NonlinearSolver.hpp:261
Scalar relaxMax() const
The greatest relaxation factor (i.e. smallest factor) allowed.
Definition NonlinearSolver.hpp:249
int linearizations() const
Number of linearizations used in all calls to step().
Definition NonlinearSolver.hpp:192
int wellIterationsLastStep() const
Number of well iterations used in all calls to step().
Definition NonlinearSolver.hpp:216
int nonlinearIterations() const
Number of full nonlinear solver iterations used in all calls to step().
Definition NonlinearSolver.hpp:196
void setParameters(const SolverParameters ¶m)
Set parameters to override those given at construction time.
Definition NonlinearSolver.hpp:265
NonlinearRelaxType relaxType() const
The relaxation type (Dampen or SOR).
Definition NonlinearSolver.hpp:257
int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition NonlinearSolver.hpp:200
int nonlinearIterationsLastStep() const
Number of nonlinear solver iterations used in the last call to step().
Definition NonlinearSolver.hpp:208
const SimulatorReportSingle & failureReport() const
return the statistics if the step() method failed
Definition NonlinearSolver.hpp:188
int wellIterations() const
Number of well iterations used in all calls to step().
Definition NonlinearSolver.hpp:204
const Model & model() const
Definition NonlinearSolver.hpp:224
void detectOscillations(const std::vector< std::vector< Scalar > > &residualHistory, const int it, bool &oscillate, bool &stagnate) const
Detect oscillation or stagnation in a given residual history.
Definition NonlinearSolver.hpp:232
NonlinearSolver(const SolverParameters ¶m, std::unique_ptr< PhysicalModel > model)
Construct solver for a given model.
Definition NonlinearSolver.hpp:111
Scalar relaxIncrement() const
The step-change size for the relaxation factor.
Definition NonlinearSolver.hpp:253
PhysicalModel & model()
Mutable reference to physical model.
Definition NonlinearSolver.hpp:228
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].
TimeStepControlInterface.
Definition TimeStepControlInterface.hpp:51
virtual bool timeStepAccepted(const double error, const double timeStepJustCompleted) const =0
For the general 3rd order controller, the internal shifting of errors and time steps happens here,...
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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
The Opm property system, traits with inheritance.
Definition NonlinearSolver.hpp:79
Definition NonlinearSolver.hpp:45
Definition NonlinearSolver.hpp:47
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34