20#ifndef OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
21#define OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
24#ifndef OPM_ADAPTIVE_TIME_STEPPING_HPP
26#include <opm/simulators/timestepping/AdaptiveTimeStepping.hpp>
27#include <opm/simulators/timestepping/AdaptiveSimulatorTimer.hpp>
30#include <dune/istl/istlexception.hh>
32#include <opm/common/Exceptions.hpp>
33#include <opm/common/ErrorMacros.hpp>
34#include <opm/common/OpmLog/OpmLog.hpp>
36#include <opm/grid/utility/StopWatch.hpp>
38#include <opm/input/eclipse/Schedule/Tuning.hpp>
40#include <opm/input/eclipse/Units/Units.hpp>
41#include <opm/input/eclipse/Units/UnitSystem.hpp>
45#include <opm/simulators/timestepping/EclTimeSteppingParams.hpp>
53#include <boost/date_time/posix_time/posix_time.hpp>
54#include <fmt/format.h>
55#include <fmt/ranges.h>
63template<
class TypeTag>
70 : time_step_control_{}
71 , restart_factor_{Parameters::Get<Parameters::SolverRestartFactor<Scalar>>()}
72 , growth_factor_{Parameters::Get<Parameters::SolverGrowthFactor<Scalar>>()}
73 , max_growth_{Parameters::Get<Parameters::SolverMaxGrowth<Scalar>>()}
75 Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60}
78 Parameters::Get<Parameters::SolverMinTimeStep<Scalar>>())}
79 , ignore_convergence_failure_{
80 Parameters::Get<Parameters::SolverContinueOnConvergenceFailure>()}
81 , solver_restart_max_{Parameters::Get<Parameters::SolverMaxRestarts>()}
82 , solver_verbose_{Parameters::Get<Parameters::SolverVerbosity>() > 0 &&
terminal_output}
83 , timestep_verbose_{Parameters::Get<Parameters::TimeStepVerbosity>() > 0 &&
terminal_output}
84 , suggested_next_timestep_{
85 (
max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>()
87 , full_timestep_initially_{Parameters::Get<Parameters::FullTimeStepInitially>()}
88 , timestep_after_event_{
89 Parameters::Get<Parameters::TimeStepAfterEventInDays<Scalar>>() * 24 * 60 * 60}
90 , use_newton_iteration_{
false}
91 , min_time_step_before_shutting_problematic_wells_{
92 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() *
unit::
day}
104template<
class TypeTag>
112 : time_step_control_{}
118 , ignore_convergence_failure_{
true}
119 , solver_restart_max_{Parameters::Get<Parameters::SolverMaxRestarts>()}
120 , solver_verbose_{Parameters::Get<Parameters::SolverVerbosity>() > 0 &&
terminal_output}
121 , timestep_verbose_{Parameters::Get<Parameters::TimeStepVerbosity>() > 0 &&
terminal_output}
122 , suggested_next_timestep_{
123 max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>() * 24 * 60 * 60
125 , full_timestep_initially_{Parameters::Get<Parameters::FullTimeStepInitially>()}
127 , use_newton_iteration_{
false}
128 , min_time_step_before_shutting_problematic_wells_{
129 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() *
unit::
day}
135template<
class TypeTag>
147 switch (this->time_step_control_type_) {
148 case TimeStepControlType::HardCodedTimeStep:
151 case TimeStepControlType::PIDAndIterationCount:
154 case TimeStepControlType::SimpleIterationCount:
157 case TimeStepControlType::PID:
160 case TimeStepControlType::General3rdOrder:
177 this->min_time_step_before_shutting_problematic_wells_ ==
181template<
class TypeTag>
183AdaptiveTimeStepping<TypeTag>::
187 detail::registerAdaptiveParameters();
198template<
class TypeTag>
199template <
class Solver>
213template<
class TypeTag>
214template<
class Serializer>
220 switch (this->time_step_control_type_) {
221 case TimeStepControlType::HardCodedTimeStep:
224 case TimeStepControlType::PIDAndIterationCount:
227 case TimeStepControlType::SimpleIterationCount:
230 case TimeStepControlType::PID:
233 case TimeStepControlType::General3rdOrder:
242 serializer(this->ignore_convergence_failure_);
250 serializer(this->min_time_step_before_shutting_problematic_wells_);
253template<
class TypeTag>
255AdaptiveTimeStepping<TypeTag>::
261template<
class TypeTag>
263AdaptiveTimeStepping<TypeTag>::
264serializationTestObjectHardcoded()
269template<
class TypeTag>
271AdaptiveTimeStepping<TypeTag>::
272serializationTestObjectPID()
277template<
class TypeTag>
279AdaptiveTimeStepping<TypeTag>::
280serializationTestObjectPIDIt()
285template<
class TypeTag>
287AdaptiveTimeStepping<TypeTag>::
288serializationTestObjectSimple()
293template<
class TypeTag>
295AdaptiveTimeStepping<TypeTag>::
296serializationTestObject3rdOrder()
302template<
class TypeTag>
304AdaptiveTimeStepping<TypeTag>::
305setSuggestedNextStep(
const double x)
307 this->suggested_next_timestep_ = x;
310template<
class TypeTag>
312AdaptiveTimeStepping<TypeTag>::
313suggestedNextStep()
const
315 return this->suggested_next_timestep_;
318template<
class TypeTag>
319const TimeStepControlInterface&
320AdaptiveTimeStepping<TypeTag>::
321timeStepControl()
const
323 return *this->time_step_control_;
327template<
class TypeTag>
329AdaptiveTimeStepping<TypeTag>::
339template<
class TypeTag>
341AdaptiveTimeStepping<TypeTag>::
344 this->restart_factor_ =
tuning.TSFCNV;
345 this->growth_factor_ =
tuning.TFDIFF;
346 this->max_growth_ =
tuning.TSFMAX;
347 this->max_time_step_ =
tuning.TSMAXZ;
349 this->timestep_after_event_ =
tuning.TMAXWC;
356template<
class TypeTag>
357template<
class T,
class Serializer>
359AdaptiveTimeStepping<TypeTag>::
363 this->time_step_control_ = std::make_unique<T>();
365 serializer(*
static_cast<T*
>(this->time_step_control_.get()));
368template<
class TypeTag>
371AdaptiveTimeStepping<TypeTag>::
374 const T* lhs =
static_cast<const T*
>(this->time_step_control_.get());
375 const T* rhs =
static_cast<const T*
>(
Rhs.time_step_control_.get());
379template<
class TypeTag>
381AdaptiveTimeStepping<TypeTag>::
386 if (this->suggested_next_timestep_ < 0) {
390 if (this->full_timestep_initially_) {
395 if (
is_event && this->timestep_after_event_ > 0) {
396 this->suggested_next_timestep_ = this->timestep_after_event_;
400template<
class TypeTag>
401template<
class Controller>
403AdaptiveTimeStepping<TypeTag>::
404serializationTestObject_()
409 result.growth_factor_ = 2.0;
411 result.max_time_step_ = 4.0;
412 result.min_time_step_ = 5.0;
413 result.ignore_convergence_failure_ =
true;
414 result.solver_restart_max_ = 6;
415 result.solver_verbose_ =
true;
416 result.timestep_verbose_ =
true;
417 result.suggested_next_timestep_ = 7.0;
418 result.full_timestep_initially_ =
true;
419 result.use_newton_iteration_ =
true;
420 result.min_time_step_before_shutting_problematic_wells_ = 9.0;
421 result.time_step_control_type_ = Controller::Type;
422 result.time_step_control_ =
423 std::make_unique<Controller>(Controller::serializationTestObject());
432template<
class TypeTag>
433void AdaptiveTimeStepping<TypeTag>::
436 std::tie(time_step_control_type_,
438 use_newton_iteration_) = detail::createController(
unitSystem);
440 if (this->growth_factor_ < 1.0) {
442 "Growth factor cannot be less than 1.");
452template<
class TypeTag>
453template<
class Solver>
454AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
468template<
class TypeTag>
469template<
class Solver>
471AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
472getAdaptiveTimerStepper()
474 return adaptive_time_stepping_;
477template<
class TypeTag>
478template<
class Solver>
480AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
483#ifdef RESERVOIR_COUPLING_ENABLED
491 return runStepOriginal_();
494 return runStepOriginal_();
502#ifdef RESERVOIR_COUPLING_ENABLED
503template<
class TypeTag>
504template<
class Solver>
506AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
507isReservoirCouplingMaster_()
const
509 return this->solver_.model().simulator().reservoirCouplingMaster() !=
nullptr;
512template<
class TypeTag>
513template<
class Solver>
515AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
516isReservoirCouplingSlave_()
const
518 return this->solver_.model().simulator().reservoirCouplingSlave() !=
nullptr;
522template<
class TypeTag>
523template<
class Solver>
525AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
528 this->adaptive_time_stepping_.maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
536template<
class TypeTag>
537template<
class Solver>
539AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
545template<
class TypeTag>
546template<
class Solver>
548AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
551 return this->adaptive_time_stepping_.max_time_step_;
554template <
class TypeTag>
555template <
class Solver>
557AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
560 const auto elapsed = this->simulator_timer_.simulationTimeElapsed();
562 const auto report_step = this->simulator_timer_.reportStepNum();
567 this->simulator_timer_.startDateTime(),
570 suggestedNextTimestep_(),
578#ifdef RESERVOIR_COUPLING_ENABLED
579template <
class TypeTag>
580template <
class Solver>
581ReservoirCouplingMaster&
582AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
583reservoirCouplingMaster_()
585 return *(this->solver_.model().simulator().reservoirCouplingMaster());
589#ifdef RESERVOIR_COUPLING_ENABLED
590template <
class TypeTag>
591template <
class Solver>
592ReservoirCouplingSlave&
593AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
594reservoirCouplingSlave_()
596 return *(this->solver_.model().simulator().reservoirCouplingSlave());
600#ifdef RESERVOIR_COUPLING_ENABLED
631template <
class TypeTag>
632template <
class Solver>
634AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
635runStepReservoirCouplingMaster_()
639 double current_time{this->simulator_timer_.simulationTimeElapsed()};
642 SimulatorReport report;
657 this->simulator_timer_.startDateTime(),
660 suggestedNextTimestep_(),
661 this->simulator_timer_.reportStepNum(),
680#ifdef RESERVOIR_COUPLING_ENABLED
681template <
class TypeTag>
682template <
class Solver>
684AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
685runStepReservoirCouplingSlave_()
689 double current_time{this->simulator_timer_.simulationTimeElapsed()};
691 SimulatorReport report;
697 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
timestep);
700 this->simulator_timer_.startDateTime(),
703 suggestedNextTimestep_(),
704 this->simulator_timer_.reportStepNum(),
724template <
class TypeTag>
725template <
class Solver>
727AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
728suggestedNextTimestep_()
const
730 return this->adaptive_time_stepping_.suggestedNextStep();
739template<
class TypeTag>
740template<
class Solver>
741AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
750 , adaptive_time_stepping_{
substepper.getAdaptiveTimerStepper()}
754template <
class TypeTag>
755template <
class Solver>
757AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
760 auto& simulator = solver_().model().simulator();
761 auto& problem = simulator.problem();
764 SimulatorReport report;
767 while (!this->substep_timer_.done()) {
771 maybeUpdateTuningAndTimeStep_();
773 const double dt = this->substep_timer_.currentStepLength();
774 if (timeStepVerbose_()) {
775 detail::logTimer(this->substep_timer_);
782 auto&
full_report = adaptive_time_stepping_.report();
789 ++this->substep_timer_;
792 auto dt_estimate = timeStepControlComputeEstimate_(
793 dt, iterations, this->substep_timer_);
800 if (this->final_step_ && this->substep_timer_.done()) {
805 report.success.output_write_time += writeOutput_();
811 report.success.converged = this->substep_timer_.done();
812 this->substep_timer_.setLastStepFailed(
false);
815 this->substep_timer_.setLastStepFailed(
true);
816 checkTimeStepMaxRestartLimit_(
restarts);
820 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>();
821 const double safetyFactor = Parameters::Get<Parameters::TimeStepControlSafetyFactor>();
841 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
843 updateSuggestedNextStep_();
853template<
class TypeTag>
854template<
class Solver>
856AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
857checkContinueOnUnconvergedSolution_(
double dt)
const
862 const auto msg = fmt::format(
863 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
864 "which is the minimum threshold given by option --solver-min-time-step\n",
867 OpmLog::problem(
msg);
872template<
class TypeTag>
873template<
class Solver>
875AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
876checkTimeStepMaxRestartLimit_(
const int restarts)
const
880 if (
restarts >= solverRestartMax_()) {
881 const auto msg = fmt::format(
882 "Solver failed to converge after cutting timestep {} times.",
restarts
884 if (solverVerbose_()) {
892template<
class TypeTag>
893template<
class Solver>
895AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
898 using Meas = UnitSystem::measure;
902 auto msg = fmt::format(
"Solver failed to converge after cutting timestep to ");
903 if (Parameters::Get<Parameters::EnableTuning>()) {
904 const UnitSystem&
unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
906 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
913 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
914 minTimeStep_() / 86400.0
917 if (solverVerbose_()) {
925template<
class TypeTag>
926template<
class Solver>
928AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
932 if (solverVerbose_()) {
933 const auto msg = fmt::format(
"{}\nTimestep chopped to {} days\n",
934 this->cause_of_failure_,
935 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
936 OpmLog::problem(
msg);
940template<
class TypeTag>
941template<
class Solver>
943AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
954 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
966 solver_().model().wellModel().forceShutWellByName(well,
967 this->substep_timer_.simulationTimeElapsed(),
977 solver_().model().wellModel().forceShutWellByName(well,
978 this->substep_timer_.simulationTimeElapsed(),
990 if (solverVerbose_()) {
991 const std::string
msg =
992 fmt::format(
"\nProblematic well(s) were shut: {}"
993 "(retrying timestep)\n",
995 OpmLog::problem(
msg);
1002template<
class TypeTag>
1003template<
class Solver>
1004boost::posix_time::ptime
1005AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1006currentDateTime_()
const
1008 return simulatorTimer_().currentDateTime();
1011template<
class TypeTag>
1012template<
class Solver>
1014AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1015getNumIterations_(
const SimulatorReportSingle &
substep_report)
const
1017 if (useNewtonIteration_()) {
1025template<
class TypeTag>
1026template<
class Solver>
1028AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1029growthFactor_()
const
1031 return this->adaptive_time_stepping_.growth_factor_;
1034template<
class TypeTag>
1035template<
class Solver>
1037AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1038ignoreConvergenceFailure_()
const
1040 return adaptive_time_stepping_.ignore_convergence_failure_;
1043template<
class TypeTag>
1044template<
class Solver>
1046AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1049 return this->adaptive_time_stepping_.max_growth_;
1052template<
class TypeTag>
1053template<
class Solver>
1055AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1058 if (timeStepVerbose_()) {
1059 std::ostringstream
ss;
1061 OpmLog::info(
ss.str());
1065template<
class TypeTag>
1066template<
class Solver>
1068AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1085template<
class TypeTag>
1086template<
class Solver>
1088AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1089maybeUpdateTuningAndTimeStep_()
1097 const auto old_value = suggestedNextTimestep_();
1098 if (this->substepper_.maybeUpdateTuning_(
this->substep_timer_.simulationTimeElapsed(),
1099 this->substep_timer_.currentStepLength(),
1100 this->substep_timer_.currentStepNum()))
1107 setTimeStep_(suggestedNextTimestep_());
1112template<
class TypeTag>
1113template<
class Solver>
1115AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1116minTimeStepBeforeClosingWells_()
const
1118 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1121template<
class TypeTag>
1122template<
class Solver>
1124AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1127 return this->adaptive_time_stepping_.min_time_step_;
1130template<
class TypeTag>
1131template<
class Solver>
1133AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1134restartFactor_()
const
1136 return this->adaptive_time_stepping_.restart_factor_;
1139template<
class TypeTag>
1140template<
class Solver>
1141SimulatorReportSingle
1142AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1153 OpmLog::debug(std::string(
"Caught Exception: ") +
e.what());
1158 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1159 if (solverVerbose_()) {
1161 OpmLog::debug(
"Overall linear iterations used: "
1166 handleFailure(
"Solver convergence failure - Iteration limit reached",
e);
1171 catch (
const ConvergenceMonitorFailure&
e) {
1178 handleFailure(
"Solver convergence failure - Numerical problem encountered",
e);
1180 catch (
const std::runtime_error&
e) {
1183 catch (
const Dune::ISTLError&
e) {
1186 catch (
const Dune::MatrixBlockError&
e) {
1193template<
class TypeTag>
1194template<
class Solver>
1196AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1199 this->substep_timer_.provideTimeStepEstimate(
dt_estimate);
1202template<
class TypeTag>
1203template<
class Solver>
1205AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1208 return this->substepper_.solver_;
1212template<
class TypeTag>
1213template<
class Solver>
1215AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1216solverRestartMax_()
const
1218 return this->adaptive_time_stepping_.solver_restart_max_;
1221template<
class TypeTag>
1222template<
class Solver>
1224AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1225setSuggestedNextStep_(
double step)
1227 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1230template <
class TypeTag>
1231template <
class Solver>
1232const SimulatorTimer&
1233AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1234simulatorTimer_()
const
1236 return this->substepper_.simulator_timer_;
1239template <
class TypeTag>
1240template <
class Solver>
1242AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1243solverVerbose_()
const
1245 return this->adaptive_time_stepping_.solver_verbose_;
1248template<
class TypeTag>
1249template<
class Solver>
1250boost::posix_time::ptime
1251AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1252startDateTime_()
const
1254 return simulatorTimer_().startDateTime();
1257template <
class TypeTag>
1258template <
class Solver>
1260AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1261suggestedNextTimestep_()
const
1263 return this->adaptive_time_stepping_.suggestedNextStep();
1266template <
class TypeTag>
1267template <
class Solver>
1269AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1270timeStepControlComputeEstimate_(
const double dt,
const int iterations,
1275 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1279template <
class TypeTag>
1280template <
class Solver>
1282AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1283timeStepVerbose_()
const
1285 return this->adaptive_time_stepping_.timestep_verbose_;
1295template <
class TypeTag>
1296template <
class Solver>
1298AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1299updateSuggestedNextStep_()
1305 if (timeStepVerbose_()) {
1306 std::ostringstream
ss;
1307 this->substep_timer_.report(
ss);
1308 ss <<
"Suggested next step size = "
1310 OpmLog::debug(
ss.str());
1315template <
class TypeTag>
1316template <
class Solver>
1318AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1319useNewtonIteration_()
const
1321 return this->adaptive_time_stepping_.use_newton_iteration_;
1324template <
class TypeTag>
1325template <
class Solver>
1327AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1332 auto& problem = solver_().model().simulator().problem();
1333 problem.writeOutput(
true);
1341template<
class TypeTag>
1342template<
class Solver>
1343AdaptiveTimeStepping<TypeTag>::
1344SolutionTimeErrorSolverWrapper<Solver>::
1345SolutionTimeErrorSolverWrapper(
const Solver& solver)
1349template<
class TypeTag>
1350template<
class Solver>
1351double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange()
const
1354 return solver_.model().relativeChange();
Definition AdaptiveTimeStepping.hpp:78
double max_growth_
factor that limits the maximum growth of a time step
Definition AdaptiveTimeStepping.hpp:252
double max_time_step_
maximal allowed time step size in days
Definition AdaptiveTimeStepping.hpp:253
bool solver_verbose_
solver verbosity
Definition AdaptiveTimeStepping.hpp:257
int solver_restart_max_
how many restart of solver are allowed
Definition AdaptiveTimeStepping.hpp:256
double timestep_after_event_
suggested size of timestep after an event
Definition AdaptiveTimeStepping.hpp:261
bool ignore_convergence_failure_
continue instead of stop when minimum time step is reached
Definition AdaptiveTimeStepping.hpp:255
TimeStepControlType time_step_control_type_
type of time step control object
Definition AdaptiveTimeStepping.hpp:248
bool full_timestep_initially_
beginning with the size of the time step from data file
Definition AdaptiveTimeStepping.hpp:260
SimulatorReport step(const SimulatorTimer &simulator_timer, Solver &solver, const bool is_event, const TuningUpdateCallback &tuning_updater)
step method that acts like the solver::step method in a sub cycle of time steps
Definition AdaptiveTimeStepping_impl.hpp:202
double growth_factor_
factor to multiply time step when solver recovered from failed convergence
Definition AdaptiveTimeStepping.hpp:251
double restart_factor_
factor to multiply time step with when solver fails to converge
Definition AdaptiveTimeStepping.hpp:250
double min_time_step_
minimal allowed time step size before throwing
Definition AdaptiveTimeStepping.hpp:254
TimeStepController time_step_control_
time step control object
Definition AdaptiveTimeStepping.hpp:249
double min_time_step_before_shutting_problematic_wells_
< shut problematic wells when time step size in days are less than this
Definition AdaptiveTimeStepping.hpp:265
bool use_newton_iteration_
use newton iteration count for adaptive time step control
Definition AdaptiveTimeStepping.hpp:262
Definition SimulatorTimer.hpp:39
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
This file provides the infrastructure to retrieve run-time parameters.
static bool compare_gt_or_eq(double a, double b)
Determines if a is greater than b within the specified tolerance.
Definition ReservoirCoupling.cpp:120
Definition SimulatorReport.hpp:122