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>
38 namespace Opm::Properties {
44 template<
class TypeTag,
class MyTypeTag>
46 using type = UndefinedProperty;
52 template<
class TypeTag,
class MyTypeTag>
54 using type = UndefinedProperty;
56 template<
class TypeTag,
class MyTypeTag>
58 using type = UndefinedProperty;
61 template<
class TypeTag>
63 using type = GetPropType<TypeTag, Scalar>;
64 static constexpr type value = 0.5;
66 template<
class TypeTag>
67 struct NewtonMaxIterations<TypeTag, TTag::FlowNonLinearSolver> {
68 static constexpr
int value = 20;
70 template<
class TypeTag>
72 static constexpr
int value = 1;
74 template<
class TypeTag>
76 static constexpr
auto value =
"dampen";
87 template <
class TypeTag,
class PhysicalModel>
90 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
102 RelaxType relaxType_;
104 double relaxIncrement_;
115 relaxMax_ = EWOMS_GET_PARAM(TypeTag, Scalar, NewtonMaxRelax);
116 maxIter_ = EWOMS_GET_PARAM(TypeTag,
int, NewtonMaxIterations);
117 minIter_ = EWOMS_GET_PARAM(TypeTag,
int, NewtonMinIterations);
119 const auto& relaxationTypeString = EWOMS_GET_PARAM(TypeTag, std::string, NewtonRelaxationType);
120 if (relaxationTypeString ==
"dampen") {
122 }
else if (relaxationTypeString ==
"sor") {
125 OPM_THROW(std::runtime_error,
"Unknown Relaxtion Type " << relaxationTypeString);
129 static void registerParameters()
131 EWOMS_REGISTER_PARAM(TypeTag, Scalar, NewtonMaxRelax,
"The maximum relaxation factor of a Newton iteration");
132 EWOMS_REGISTER_PARAM(TypeTag,
int, NewtonMaxIterations,
"The maximum number of Newton iterations per time step");
133 EWOMS_REGISTER_PARAM(TypeTag,
int, NewtonMinIterations,
"The minimum number of Newton iterations per time step");
134 EWOMS_REGISTER_PARAM(TypeTag, std::string, NewtonRelaxationType,
"The type of relaxation used by Newton method");
142 relaxIncrement_ = 0.1;
163 std::unique_ptr<PhysicalModel>
model)
165 , model_(std::move(
model))
167 , nonlinearIterations_(0)
168 , linearIterations_(0)
170 , nonlinearIterationsLast_(0)
171 , linearIterationsLast_(0)
172 , wellIterationsLast_(0)
175 OPM_THROW(std::logic_error,
"Must provide a non-null model argument for NonlinearSolver.");
187 report += model_->prepareStep(timer);
194 bool converged =
false;
202 auto iterReport = model_->nonlinearIteration(iteration, timer, *
this);
204 report += iterReport;
205 report.converged = iterReport.converged;
207 converged = report.converged;
213 failureReport_ = report;
214 failureReport_ += model_->failureReport();
218 while ( (!converged && (iteration <=
maxIter())) || (iteration <=
minIter()));
221 failureReport_ = report;
223 std::string msg =
"Solver convergence failure - Failed to complete a time step within " + std::to_string(
maxIter()) +
" iterations.";
224 OPM_THROW_NOLOG(TooManyIterations, msg);
228 report += model_->afterStep(timer);
229 report.converged =
true;
235 {
return failureReport_; }
239 {
return linearizations_; }
243 {
return nonlinearIterations_; }
247 {
return linearIterations_; }
251 {
return wellIterations_; }
255 {
return nonlinearIterationsLast_; }
259 {
return linearIterationsLast_; }
263 {
return wellIterationsLast_; }
265 std::vector<std::vector<double> >
266 computeFluidInPlace(
const std::vector<int>& fipnum)
const
267 {
return model_->computeFluidInPlace(fipnum); }
279 const int it,
bool& oscillate,
bool& stagnate)
const
293 int oscillatePhase = 0;
294 const std::vector<double>& F0 = residualHistory[it];
295 const std::vector<double>& F1 = residualHistory[it - 1];
296 const std::vector<double>& F2 = residualHistory[it - 2];
297 for (
int p= 0; p < model_->numPhases(); ++p){
298 const double d1 = std::abs((F0[p] - F2[p]) / F0[p]);
299 const double d2 = std::abs((F0[p] - F1[p]) / F0[p]);
305 stagnate = (stagnate && !(std::abs((F1[p] - F2[p]) / F2[p]) > 1.0e-3));
308 oscillate = (oscillatePhase > 1);
314 template <
class BVector>
320 BVector tempDxOld = dxOld;
329 for (i = 0; i < dx.size(); ++i) {
339 for (i = 0; i < dx.size(); ++i) {
341 tempDxOld[i] *= (1.-omega);
342 dx[i] += tempDxOld[i];
347 OPM_THROW(std::runtime_error,
"Can only handle Dampen and SOR relaxation type.");
355 {
return param_.relaxMax_; }
359 {
return param_.relaxIncrement_; }
363 {
return param_.relaxType_; }
367 {
return param_.relaxRelTol_; }
371 {
return param_.maxIter_; }
375 {
return param_.minIter_; }
384 SolverParameters param_;
385 std::unique_ptr<PhysicalModel> model_;
387 int nonlinearIterations_;
388 int linearIterations_;
390 int nonlinearIterationsLast_;
391 int linearIterationsLast_;
392 int wellIterationsLast_;
A nonlinear solver class suitable for general fully-implicit models, as well as pressure,...
Definition: NonlinearSolverEbos.hpp:89
double relaxMax() const
The greatest relaxation factor (i.e. smallest factor) allowed.
Definition: NonlinearSolverEbos.hpp:354
int maxIter() const
The maximum number of nonlinear iterations allowed.
Definition: NonlinearSolverEbos.hpp:370
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:278
PhysicalModel & model()
Mutable reference to physical model.
Definition: NonlinearSolverEbos.hpp:274
int wellIterationsLastStep() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:262
const SimulatorReportSingle & failureReport() const
return the statistics if the step() method failed
Definition: NonlinearSolverEbos.hpp:234
double relaxIncrement() const
The step-change size for the relaxation factor.
Definition: NonlinearSolverEbos.hpp:358
int linearIterationsLastStep() const
Number of linear solver iterations used in the last call to step().
Definition: NonlinearSolverEbos.hpp:258
int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:246
enum RelaxType relaxType() const
The relaxation type (Dampen or SOR).
Definition: NonlinearSolverEbos.hpp:362
int linearizations() const
Number of linearizations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:238
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const double omega) const
Apply a stabilization to dx, depending on dxOld and relaxation parameters.
Definition: NonlinearSolverEbos.hpp:315
int minIter() const
The minimum number of nonlinear iterations allowed.
Definition: NonlinearSolverEbos.hpp:374
int nonlinearIterationsLastStep() const
Number of nonlinear solver iterations used in the last call to step().
Definition: NonlinearSolverEbos.hpp:254
double relaxRelTol() const
The relaxation relative tolerance.
Definition: NonlinearSolverEbos.hpp:366
int wellIterations() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:250
const PhysicalModel & model() const
Reference to physical model.
Definition: NonlinearSolverEbos.hpp:270
NonlinearSolverEbos(const SolverParameters ¶m, std::unique_ptr< PhysicalModel > model)
Construct solver for a given model.
Definition: NonlinearSolverEbos.hpp:162
void setParameters(const SolverParameters ¶m)
Set parameters to override those given at construction time.
Definition: NonlinearSolverEbos.hpp:378
int nonlinearIterations() const
Number of full nonlinear solver iterations used in all calls to step().
Definition: NonlinearSolverEbos.hpp:242
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:56
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27
Definition: NonlinearSolverEbos.hpp:101
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:31