20#ifndef OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
21#define OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
24#ifndef OPM_ADAPTIVE_TIME_STEPPING_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>
53#include <boost/date_time/posix_time/posix_time.hpp>
54#include <fmt/format.h>
55#include <fmt/ranges.h>
63template<
class TypeTag>
67 const double max_next_tstep,
68 const bool terminal_output
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}
77 unit_system.to_si(UnitSystem::measure::time,
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>()
86 : max_next_tstep) * 24 * 60 * 60}
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}
100template<
class TypeTag>
103 const Tuning& tuning,
104 const UnitSystem& unit_system,
106 const bool terminal_output
108 : time_step_control_{}
109 , restart_factor_{tuning.TSFCNV}
110 , growth_factor_{tuning.TFDIFF}
111 , max_growth_{tuning.TSFMAX}
112 , max_time_step_{tuning.TSMAXZ}
113 , min_time_step_{tuning.TSMINZ}
114 , ignore_convergence_failure_{true}
115 , solver_restart_max_{Parameters::
Get<Parameters::SolverMaxRestarts>()}
116 , solver_verbose_{Parameters::
Get<Parameters::SolverVerbosity>() > 0 && terminal_output}
117 , timestep_verbose_{Parameters::
Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output}
118 , suggested_next_timestep_{
119 max_next_tstep <= 0 ? Parameters::
Get<Parameters::InitialTimeStepInDays>() * 24 * 60 * 60
121 , full_timestep_initially_{Parameters::
Get<Parameters::FullTimeStepInitially>()}
122 , timestep_after_event_{tuning.TMAXWC}
123 , use_newton_iteration_{false}
124 , min_time_step_before_shutting_problematic_wells_{
125 Parameters::
Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
131template<
class TypeTag>
143 switch (this->time_step_control_type_) {
145 result = castAndComp<HardcodedTimeStepControl>(rhs);
148 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
151 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
154 result = castAndComp<PIDTimeStepControl>(rhs);
157 result = castAndComp<General3rdOrderController>(rhs);
173 this->min_time_step_before_shutting_problematic_wells_ ==
177template<
class TypeTag>
182 registerEclTimeSteppingParameters<Scalar>();
191template<
class TypeTag>
192template <
class Solver>
200 SubStepper<Solver> sub_stepper{
201 *
this, simulator_timer, solver, is_event, tuning_updater,
203 return sub_stepper.run();
206template<
class TypeTag>
207template<
class Serializer>
212 serializer(this->time_step_control_type_);
213 switch (this->time_step_control_type_) {
215 allocAndSerialize<HardcodedTimeStepControl>(serializer);
218 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
221 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
224 allocAndSerialize<PIDTimeStepControl>(serializer);
227 allocAndSerialize<General3rdOrderController>(serializer);
230 serializer(this->restart_factor_);
231 serializer(this->growth_factor_);
232 serializer(this->max_growth_);
233 serializer(this->max_time_step_);
234 serializer(this->min_time_step_);
235 serializer(this->ignore_convergence_failure_);
236 serializer(this->solver_restart_max_);
237 serializer(this->solver_verbose_);
238 serializer(this->timestep_verbose_);
239 serializer(this->suggested_next_timestep_);
240 serializer(this->full_timestep_initially_);
241 serializer(this->timestep_after_event_);
242 serializer(this->use_newton_iteration_);
243 serializer(this->min_time_step_before_shutting_problematic_wells_);
246template<
class TypeTag>
254template<
class TypeTag>
259 return serializationTestObject_<HardcodedTimeStepControl>();
262template<
class TypeTag>
267 return serializationTestObject_<PIDTimeStepControl>();
270template<
class TypeTag>
275 return serializationTestObject_<PIDAndIterationCountTimeStepControl>();
278template<
class TypeTag>
283 return serializationTestObject_<SimpleIterationCountTimeStepControl>();
286template<
class TypeTag>
291 return serializationTestObject_<General3rdOrderController>();
295template<
class TypeTag>
300 this->suggested_next_timestep_ = x;
303template<
class TypeTag>
308 return this->suggested_next_timestep_;
311template<
class TypeTag>
316 return *this->time_step_control_;
320template<
class TypeTag>
327 if (max_next_tstep > 0) {
328 this->suggested_next_timestep_ = max_next_tstep;
332template<
class TypeTag>
335updateTUNING(
double max_next_tstep,
const Tuning& tuning)
337 this->restart_factor_ = tuning.TSFCNV;
338 this->growth_factor_ = tuning.TFDIFF;
339 this->max_growth_ = tuning.TSFMAX;
340 this->max_time_step_ = tuning.TSMAXZ;
341 updateNEXTSTEP(max_next_tstep);
342 this->timestep_after_event_ = tuning.TMAXWC;
349template<
class TypeTag>
350template<
class T,
class Serializer>
355 if (!serializer.isSerializing()) {
356 this->time_step_control_ = std::make_unique<T>();
358 serializer(*
static_cast<T*
>(this->time_step_control_.get()));
361template<
class TypeTag>
364AdaptiveTimeStepping<TypeTag>::
365castAndComp(
const AdaptiveTimeStepping<TypeTag>& Rhs)
const
367 const T* lhs =
static_cast<const T*
>(this->time_step_control_.get());
368 const T* rhs =
static_cast<const T*
>(Rhs.time_step_control_.get());
372template<
class TypeTag>
374AdaptiveTimeStepping<TypeTag>::
375maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
const double original_time_step,
379 if (this->suggested_next_timestep_ < 0) {
380 this->suggested_next_timestep_ = this->restart_factor_ * original_time_step;
383 if (this->full_timestep_initially_) {
384 this->suggested_next_timestep_ = original_time_step;
388 if (is_event && this->timestep_after_event_ > 0) {
389 this->suggested_next_timestep_ = this->timestep_after_event_;
393template<
class TypeTag>
394template<
class Controller>
395AdaptiveTimeStepping<TypeTag>
396AdaptiveTimeStepping<TypeTag>::
397serializationTestObject_()
399 AdaptiveTimeStepping<TypeTag> result;
401 result.restart_factor_ = 1.0;
402 result.growth_factor_ = 2.0;
403 result.max_growth_ = 3.0;
404 result.max_time_step_ = 4.0;
405 result.min_time_step_ = 5.0;
406 result.ignore_convergence_failure_ =
true;
407 result.solver_restart_max_ = 6;
408 result.solver_verbose_ =
true;
409 result.timestep_verbose_ =
true;
410 result.suggested_next_timestep_ = 7.0;
411 result.full_timestep_initially_ =
true;
412 result.use_newton_iteration_ =
true;
413 result.min_time_step_before_shutting_problematic_wells_ = 9.0;
414 result.time_step_control_type_ = Controller::Type;
415 result.time_step_control_ =
416 std::make_unique<Controller>(Controller::serializationTestObject());
425template<
class TypeTag>
427init_(
const UnitSystem& unitSystem)
429 std::tie(time_step_control_type_,
433 if (this->growth_factor_ < 1.0) {
434 OPM_THROW(std::runtime_error,
435 "Growth factor cannot be less than 1.");
445template<
class TypeTag>
446template<
class Solver>
452 const TuningUpdateCallback& tuning_updater)
453 : adaptive_time_stepping_{adaptive_time_stepping}
454 , simulator_timer_{simulator_timer}
456 , is_event_{is_event}
457 , tuning_updater_{tuning_updater}
461template<
class TypeTag>
462template<
class Solver>
463AdaptiveTimeStepping<TypeTag>&
464AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
465getAdaptiveTimerStepper()
467 return adaptive_time_stepping_;
470template<
class TypeTag>
471template<
class Solver>
473AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
476#ifdef RESERVOIR_COUPLING_ENABLED
477 if (isReservoirCouplingSlave_() && reservoirCouplingSlave_().activated()) {
478 return runStepReservoirCouplingSlave_();
480 else if (isReservoirCouplingMaster_() && reservoirCouplingMaster_().activated()) {
481 return runStepReservoirCouplingMaster_();
484 return runStepOriginal_();
487 return runStepOriginal_();
495#ifdef RESERVOIR_COUPLING_ENABLED
496template<
class TypeTag>
497template<
class Solver>
499AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
500isReservoirCouplingMaster_()
const
502 return this->solver_.model().simulator().reservoirCouplingMaster() !=
nullptr;
505template<
class TypeTag>
506template<
class Solver>
508AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
509isReservoirCouplingSlave_()
const
511 return this->solver_.model().simulator().reservoirCouplingSlave() !=
nullptr;
515template<
class TypeTag>
516template<
class Solver>
518AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
519maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
const double original_time_step)
521 this->adaptive_time_stepping_.maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
522 original_time_step, this->is_event_
529template<
class TypeTag>
530template<
class Solver>
532AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
533maybeUpdateTuning_(
double elapsed,
double dt,
int sub_step_number)
const
535 return this->tuning_updater_(elapsed, dt, sub_step_number);
538template<
class TypeTag>
539template<
class Solver>
541AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
544 return this->adaptive_time_stepping_.max_time_step_;
547template <
class TypeTag>
548template <
class Solver>
550AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
553 const auto elapsed = this->simulator_timer_.simulationTimeElapsed();
554 const auto original_time_step = this->simulator_timer_.currentStepLength();
555 const auto report_step = this->simulator_timer_.reportStepNum();
556 maybeUpdateTuning_(elapsed, original_time_step, report_step);
557 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(original_time_step);
559 AdaptiveSimulatorTimer substep_timer{
560 this->simulator_timer_.startDateTime(),
563 suggestedNextTimestep_(),
567 SubStepIteration<Solver> substepIteration{*
this, substep_timer, original_time_step,
true};
568 return substepIteration.run();
571#ifdef RESERVOIR_COUPLING_ENABLED
572template <
class TypeTag>
573template <
class Solver>
574ReservoirCouplingMaster&
575AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
576reservoirCouplingMaster_()
578 return *(this->solver_.model().simulator().reservoirCouplingMaster());
582#ifdef RESERVOIR_COUPLING_ENABLED
583template <
class TypeTag>
584template <
class Solver>
585ReservoirCouplingSlave&
586AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
587reservoirCouplingSlave_()
589 return *(this->solver_.model().simulator().reservoirCouplingSlave());
593#ifdef RESERVOIR_COUPLING_ENABLED
624template <
class TypeTag>
625template <
class Solver>
627AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
628runStepReservoirCouplingMaster_()
631 const double original_time_step = this->simulator_timer_.currentStepLength();
632 double current_time{this->simulator_timer_.simulationTimeElapsed()};
633 double step_end_time = current_time + original_time_step;
634 auto current_step_length = original_time_step;
635 SimulatorReport report;
637 reservoirCouplingMaster_().receiveNextReportDateFromSlaves();
638 if (iteration == 0) {
639 maybeUpdateTuning_(current_time, current_step_length, 0);
641 current_step_length = reservoirCouplingMaster_().maybeChopSubStep(
642 current_step_length, current_time);
643 reservoirCouplingMaster_().sendNextTimeStepToSlaves(current_step_length);
644 if (iteration == 0) {
645 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
647 AdaptiveSimulatorTimer substep_timer{
648 this->simulator_timer_.startDateTime(),
651 suggestedNextTimestep_(),
652 this->simulator_timer_.reportStepNum(),
656 current_time + current_step_length, step_end_time
658 SubStepIteration<Solver> substepIteration{*
this, substep_timer, current_step_length, final_step};
659 const auto sub_steps_report = substepIteration.run();
660 report += sub_steps_report;
661 current_time += current_step_length;
671#ifdef RESERVOIR_COUPLING_ENABLED
672template <
class TypeTag>
673template <
class Solver>
675AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
676runStepReservoirCouplingSlave_()
679 const double original_time_step = this->simulator_timer_.currentStepLength();
680 double current_time{this->simulator_timer_.simulationTimeElapsed()};
681 double step_end_time = current_time + original_time_step;
682 SimulatorReport report;
684 reservoirCouplingSlave_().sendNextReportDateToMasterProcess();
685 const auto timestep = reservoirCouplingSlave_().receiveNextTimeStepFromMaster();
686 if (iteration == 0) {
687 maybeUpdateTuning_(current_time, original_time_step, 0);
688 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
690 AdaptiveSimulatorTimer substep_timer{
691 this->simulator_timer_.startDateTime(),
694 suggestedNextTimestep_(),
695 this->simulator_timer_.reportStepNum(),
699 current_time + timestep, step_end_time
701 SubStepIteration<Solver> substepIteration{*
this, substep_timer, timestep, final_step};
702 const auto sub_steps_report = substepIteration.run();
703 report += sub_steps_report;
704 current_time += timestep;
715template <
class TypeTag>
716template <
class Solver>
718AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
719suggestedNextTimestep_()
const
721 return this->adaptive_time_stepping_.suggestedNextStep();
730template<
class TypeTag>
731template<
class Solver>
732AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
733SubStepIteration(SubStepper<Solver>& substepper,
734 AdaptiveSimulatorTimer& substep_timer,
735 const double original_time_step,
737 : substepper_{substepper}
738 , substep_timer_{substep_timer}
739 , original_time_step_{original_time_step}
740 , final_step_{final_step}
741 , adaptive_time_stepping_{substepper.getAdaptiveTimerStepper()}
745template <
class TypeTag>
746template <
class Solver>
748AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
751 auto& simulator = solver_().model().simulator();
752 auto& problem = simulator.problem();
755 SimulatorReport report;
758 while (!this->substep_timer_.done()) {
762 maybeUpdateTuningAndTimeStep_();
764 const double dt = this->substep_timer_.currentStepLength();
765 if (timeStepVerbose_()) {
769 const auto substep_report = runSubStep_();
772 problem.setSubStepReport(substep_report);
773 auto& full_report = adaptive_time_stepping_.report();
774 full_report += substep_report;
775 problem.setSimulationReport(full_report);
777 report += substep_report;
779 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
780 ++this->substep_timer_;
782 const int iterations = getNumIterations_(substep_report);
783 auto dt_estimate = timeStepControlComputeEstimate_(
784 dt, iterations, this->substep_timer_);
786 assert(dt_estimate > 0);
787 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
790 maybeReportSubStep_(substep_report);
791 if (this->final_step_ && this->substep_timer_.done()) {
796 report.success.output_write_time += writeOutput_();
800 setTimeStep_(dt_estimate);
802 report.success.converged = this->substep_timer_.done();
803 this->substep_timer_.setLastStepFailed(
false);
806 this->substep_timer_.setLastStepFailed(
true);
807 checkTimeStepMaxRestartLimit_(restarts);
809 double new_time_step = restartFactor_() * dt;
810 if (substep_report.time_step_rejected) {
811 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>();
812 const double safetyFactor = Parameters::Get<Parameters::TimeStepControlSafetyFactor>();
813 const double temp_time_step = std::sqrt(safetyFactor * tol / solver_().model().relativeChange()) * dt;
814 if (temp_time_step < dt) {
815 new_time_step = temp_time_step;
818 checkTimeStepMinLimit_(new_time_step);
819 bool wells_shut =
false;
820 if (new_time_step > minTimeStepBeforeClosingWells_()) {
821 chopTimeStep_(new_time_step);
823 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
832 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
834 updateSuggestedNextStep_();
844template<
class TypeTag>
845template<
class Solver>
847AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
848checkContinueOnUnconvergedSolution_(
double dt)
const
850 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
851 if (continue_on_uncoverged_solution && solverVerbose_()) {
853 const auto msg = fmt::format(
854 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
855 "which is the minimum threshold given by option --solver-min-time-step\n",
858 OpmLog::problem(msg);
860 return continue_on_uncoverged_solution;
863template<
class TypeTag>
864template<
class Solver>
866AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
867checkTimeStepMaxRestartLimit_(
const int restarts)
const
871 if (restarts >= solverRestartMax_()) {
872 const auto msg = fmt::format(
873 "Solver failed to converge after cutting timestep {} times.", restarts
875 if (solverVerbose_()) {
879 throw TimeSteppingBreakdown{msg};
883template<
class TypeTag>
884template<
class Solver>
886AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
887checkTimeStepMinLimit_(
const double new_time_step)
const
889 using Meas = UnitSystem::measure;
892 if (new_time_step < minTimeStep_()) {
893 auto msg = fmt::format(
"Solver failed to converge after cutting timestep to ");
894 if (Parameters::Get<Parameters::EnableTuning>()) {
895 const UnitSystem& unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
897 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
898 unit_system.from_si(Meas::time, minTimeStep_()),
899 unit_system.name(Meas::time)
904 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
905 minTimeStep_() / 86400.0
908 if (solverVerbose_()) {
912 throw TimeSteppingBreakdown{msg};
916template<
class TypeTag>
917template<
class Solver>
919AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
920chopTimeStep_(
const double new_time_step)
922 setTimeStep_(new_time_step);
923 if (solverVerbose_()) {
924 const auto msg = fmt::format(
"{}\nTimestep chopped to {} days\n",
925 this->cause_of_failure_,
926 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
927 OpmLog::problem(msg);
931template<
class TypeTag>
932template<
class Solver>
934AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
935chopTimeStepOrCloseFailingWells_(
const double new_time_step)
937 bool wells_shut =
false;
944 const bool requireRepeatedFailures =
945 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
946 const std::set<std::string> failing_wells =
949 if (failing_wells.empty()) {
951 chopTimeStep_(new_time_step);
954 std::vector<std::string> shut_wells;
955 for (
const auto& well : failing_wells) {
956 const bool was_shut =
957 solver_().model().wellModel().forceShutWellByName(well,
958 this->substep_timer_.simulationTimeElapsed(),
961 shut_wells.push_back(well);
965 if (shut_wells.empty()) {
966 for (
const auto& well : failing_wells) {
967 const bool was_shut =
968 solver_().model().wellModel().forceShutWellByName(well,
969 this->substep_timer_.simulationTimeElapsed(),
972 shut_wells.push_back(well);
977 if (shut_wells.empty()) {
978 chopTimeStep_(new_time_step);
981 if (solverVerbose_()) {
982 const std::string msg =
983 fmt::format(
"\nProblematic well(s) were shut: {}"
984 "(retrying timestep)\n",
985 fmt::join(shut_wells,
" "));
986 OpmLog::problem(msg);
993template<
class TypeTag>
994template<
class Solver>
995boost::posix_time::ptime
996AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
997currentDateTime_()
const
999 return simulatorTimer_().currentDateTime();
1002template<
class TypeTag>
1003template<
class Solver>
1005AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1006getNumIterations_(
const SimulatorReportSingle &substep_report)
const
1008 if (useNewtonIteration_()) {
1009 return substep_report.total_newton_iterations;
1012 return substep_report.total_linear_iterations;
1016template<
class TypeTag>
1017template<
class Solver>
1019AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1020growthFactor_()
const
1022 return this->adaptive_time_stepping_.growth_factor_;
1025template<
class TypeTag>
1026template<
class Solver>
1028AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1029ignoreConvergenceFailure_()
const
1031 return adaptive_time_stepping_.ignore_convergence_failure_;
1034template<
class TypeTag>
1035template<
class Solver>
1037AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1040 return this->adaptive_time_stepping_.max_growth_;
1043template<
class TypeTag>
1044template<
class Solver>
1046AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1047maybeReportSubStep_(SimulatorReportSingle substep_report)
const
1049 if (timeStepVerbose_()) {
1050 std::ostringstream ss;
1051 substep_report.reportStep(ss);
1052 OpmLog::info(ss.str());
1056template<
class TypeTag>
1057template<
class Solver>
1059AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1060maybeRestrictTimeStepGrowth_(
const double dt,
double dt_estimate,
const int restarts)
const
1063 dt_estimate = std::min(dt_estimate,
double(maxGrowth_() * dt));
1064 assert(dt_estimate > 0);
1067 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1076template<
class TypeTag>
1077template<
class Solver>
1079AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1080maybeUpdateTuningAndTimeStep_()
1088 const auto old_value = suggestedNextTimestep_();
1089 if (this->substepper_.maybeUpdateTuning_(this->substep_timer_.simulationTimeElapsed(),
1090 this->substep_timer_.currentStepLength(),
1091 this->substep_timer_.currentStepNum()))
1098 setTimeStep_(suggestedNextTimestep_());
1099 setSuggestedNextStep_(old_value);
1103template<
class TypeTag>
1104template<
class Solver>
1106AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1107minTimeStepBeforeClosingWells_()
const
1109 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1112template<
class TypeTag>
1113template<
class Solver>
1115AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1118 return this->adaptive_time_stepping_.min_time_step_;
1121template<
class TypeTag>
1122template<
class Solver>
1124AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1125restartFactor_()
const
1127 return this->adaptive_time_stepping_.restart_factor_;
1130template<
class TypeTag>
1131template<
class Solver>
1132SimulatorReportSingle
1133AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1136 SimulatorReportSingle substep_report;
1138 auto handleFailure = [
this, &substep_report]
1139 (
const std::string& failure_reason,
const std::exception& e,
bool log_exception =
true)
1141 substep_report = solver_().failureReport();
1142 this->cause_of_failure_ = failure_reason;
1143 if (log_exception && solverVerbose_()) {
1144 OpmLog::debug(std::string(
"Caught Exception: ") + e.what());
1149 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1150 if (solverVerbose_()) {
1152 OpmLog::debug(
"Overall linear iterations used: "
1156 catch (
const TooManyIterations& e) {
1157 handleFailure(
"Solver convergence failure - Iteration limit reached", e);
1159 catch (
const TimeSteppingBreakdown& e) {
1160 handleFailure(e.what(), e);
1162 catch (
const ConvergenceMonitorFailure& e) {
1163 handleFailure(
"Convergence monitor failure", e,
false);
1165 catch (
const LinearSolverProblem& e) {
1166 handleFailure(
"Linear solver convergence failure", e);
1168 catch (
const NumericalProblem& e) {
1169 handleFailure(
"Solver convergence failure - Numerical problem encountered", e);
1171 catch (
const std::runtime_error& e) {
1172 handleFailure(
"Runtime error encountered", e);
1174 catch (
const Dune::ISTLError& e) {
1175 handleFailure(
"ISTL error - Time step too large", e);
1177 catch (
const Dune::MatrixBlockError& e) {
1178 handleFailure(
"Matrix block error", e);
1181 return substep_report;
1184template<
class TypeTag>
1185template<
class Solver>
1187AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1188setTimeStep_(
double dt_estimate)
1190 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1193template<
class TypeTag>
1194template<
class Solver>
1196AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1199 return this->substepper_.solver_;
1203template<
class TypeTag>
1204template<
class Solver>
1206AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1207solverRestartMax_()
const
1209 return this->adaptive_time_stepping_.solver_restart_max_;
1212template<
class TypeTag>
1213template<
class Solver>
1215AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1216setSuggestedNextStep_(
double step)
1218 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1221template <
class TypeTag>
1222template <
class Solver>
1223const SimulatorTimer&
1224AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1225simulatorTimer_()
const
1227 return this->substepper_.simulator_timer_;
1230template <
class TypeTag>
1231template <
class Solver>
1233AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1234solverVerbose_()
const
1236 return this->adaptive_time_stepping_.solver_verbose_;
1239template<
class TypeTag>
1240template<
class Solver>
1241boost::posix_time::ptime
1242AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1243startDateTime_()
const
1245 return simulatorTimer_().startDateTime();
1248template <
class TypeTag>
1249template <
class Solver>
1251AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1252suggestedNextTimestep_()
const
1254 return this->adaptive_time_stepping_.suggestedNextStep();
1257template <
class TypeTag>
1258template <
class Solver>
1260AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1261timeStepControlComputeEstimate_(
const double dt,
const int iterations,
1262 const AdaptiveSimulatorTimer& substepTimer)
const
1265 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1266 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1267 dt, iterations, relative_change, substepTimer);
1270template <
class TypeTag>
1271template <
class Solver>
1273AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1274timeStepVerbose_()
const
1276 return this->adaptive_time_stepping_.timestep_verbose_;
1286template <
class TypeTag>
1287template <
class Solver>
1289AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1290updateSuggestedNextStep_()
1292 auto suggested_next_step = this->substep_timer_.currentStepLength();
1293 if (! std::isfinite(suggested_next_step)) {
1294 suggested_next_step = this->original_time_step_;
1296 if (timeStepVerbose_()) {
1297 std::ostringstream ss;
1298 this->substep_timer_.report(ss);
1299 ss <<
"Suggested next step size = "
1300 << unit::convert::to(suggested_next_step, unit::day) <<
" (days)" << std::endl;
1301 OpmLog::debug(ss.str());
1303 setSuggestedNextStep_(suggested_next_step);
1306template <
class TypeTag>
1307template <
class Solver>
1309AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1310useNewtonIteration_()
const
1312 return this->adaptive_time_stepping_.use_newton_iteration_;
1315template <
class TypeTag>
1316template <
class Solver>
1318AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1321 time::StopWatch perf_timer;
1323 auto& problem = solver_().model().simulator().problem();
1324 problem.writeOutput(
true);
1325 return perf_timer.secsSinceStart();
1332template<
class TypeTag>
1333template<
class Solver>
1334AdaptiveTimeStepping<TypeTag>::
1335SolutionTimeErrorSolverWrapper<Solver>::
1336SolutionTimeErrorSolverWrapper(
const Solver& solver)
1340template<
class TypeTag>
1341template<
class Solver>
1342double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange()
const
1345 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
void init_(const UnitSystem &unitSystem)
Definition: AdaptiveTimeStepping_impl.hpp:427
void setSuggestedNextStep(const double x)
Definition: AdaptiveTimeStepping_impl.hpp:298
double suggestedNextStep() const
Definition: AdaptiveTimeStepping_impl.hpp:306
std::function< bool(const double, const double, const int)> TuningUpdateCallback
Definition: AdaptiveTimeStepping.hpp:80
static AdaptiveTimeStepping< TypeTag > serializationTestObjectSimple()
Definition: AdaptiveTimeStepping_impl.hpp:281
bool ignore_convergence_failure_
continue instead of stop when minimum time step is reached
Definition: AdaptiveTimeStepping.hpp:255
void serializeOp(Serializer &serializer)
Definition: AdaptiveTimeStepping_impl.hpp:210
void updateTUNING(double max_next_tstep, const Tuning &tuning)
Definition: AdaptiveTimeStepping_impl.hpp:335
TimeStepControlType time_step_control_type_
type of time step control object
Definition: AdaptiveTimeStepping.hpp:248
const TimeStepControlInterface & timeStepControl() const
Definition: AdaptiveTimeStepping_impl.hpp:314
bool full_timestep_initially_
beginning with the size of the time step from data file
Definition: AdaptiveTimeStepping.hpp:260
bool operator==(const AdaptiveTimeStepping< TypeTag > &rhs)
Definition: AdaptiveTimeStepping_impl.hpp:134
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:195
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
void updateNEXTSTEP(double max_next_tstep)
Definition: AdaptiveTimeStepping_impl.hpp:323
static AdaptiveTimeStepping< TypeTag > serializationTestObjectHardcoded()
Definition: AdaptiveTimeStepping_impl.hpp:257
AdaptiveTimeStepping()=default
TimeStepController time_step_control_
time step control object
Definition: AdaptiveTimeStepping.hpp:249
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPIDIt()
Definition: AdaptiveTimeStepping_impl.hpp:273
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
static AdaptiveTimeStepping< TypeTag > serializationTestObject3rdOrder()
Definition: AdaptiveTimeStepping_impl.hpp:289
SimulatorReport & report()
Definition: AdaptiveTimeStepping_impl.hpp:249
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPID()
Definition: AdaptiveTimeStepping_impl.hpp:265
static void registerParameters()
Definition: AdaptiveTimeStepping_impl.hpp:180
bool use_newton_iteration_
use newton iteration count for adaptive time step control
Definition: AdaptiveTimeStepping.hpp:262
Definition: SimulatorTimer.hpp:39
Definition: TimeStepControlInterface.hpp:51
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:185
void logTimer(const AdaptiveSimulatorTimer &substep_timer)
void registerAdaptiveParameters()
std::set< std::string > consistentlyFailingWells(const std::vector< StepReport > &sr, bool requireRepeatedFailures)
std::tuple< TimeStepControlType, std::unique_ptr< TimeStepControlInterface >, bool > createController(const UnitSystem &unitSystem)
Definition: blackoilboundaryratevector.hh:39
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
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: SimulatorReport.hpp:122