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}
104template<
class TypeTag>
107 const Tuning& tuning,
108 const UnitSystem& unit_system,
110 const bool terminal_output
112 : time_step_control_{}
113 , restart_factor_{tuning.TSFCNV}
114 , growth_factor_{tuning.TFDIFF}
115 , max_growth_{tuning.TSFMAX}
116 , max_time_step_{tuning.TSMAXZ}
117 , min_time_step_{tuning.TSMINZ}
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>()}
126 , timestep_after_event_{tuning.TMAXWC}
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_) {
149 result = castAndComp<HardcodedTimeStepControl>(rhs);
152 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
155 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
158 result = castAndComp<PIDTimeStepControl>(rhs);
161 result = castAndComp<General3rdOrderController>(rhs);
177 this->min_time_step_before_shutting_problematic_wells_ ==
181template<
class TypeTag>
186 registerEclTimeSteppingParameters<Scalar>();
198template<
class TypeTag>
199template <
class Solver>
207 SubStepper<Solver> sub_stepper{
208 *
this, simulator_timer, solver, is_event, tuning_updater,
210 return sub_stepper.run();
213template<
class TypeTag>
214template<
class Serializer>
219 serializer(this->time_step_control_type_);
220 switch (this->time_step_control_type_) {
222 allocAndSerialize<HardcodedTimeStepControl>(serializer);
225 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
228 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
231 allocAndSerialize<PIDTimeStepControl>(serializer);
234 allocAndSerialize<General3rdOrderController>(serializer);
237 serializer(this->restart_factor_);
238 serializer(this->growth_factor_);
239 serializer(this->max_growth_);
240 serializer(this->max_time_step_);
241 serializer(this->min_time_step_);
242 serializer(this->ignore_convergence_failure_);
243 serializer(this->solver_restart_max_);
244 serializer(this->solver_verbose_);
245 serializer(this->timestep_verbose_);
246 serializer(this->suggested_next_timestep_);
247 serializer(this->full_timestep_initially_);
248 serializer(this->timestep_after_event_);
249 serializer(this->use_newton_iteration_);
250 serializer(this->min_time_step_before_shutting_problematic_wells_);
253template<
class TypeTag>
261template<
class TypeTag>
266 return serializationTestObject_<HardcodedTimeStepControl>();
269template<
class TypeTag>
274 return serializationTestObject_<PIDTimeStepControl>();
277template<
class TypeTag>
282 return serializationTestObject_<PIDAndIterationCountTimeStepControl>();
285template<
class TypeTag>
290 return serializationTestObject_<SimpleIterationCountTimeStepControl>();
293template<
class TypeTag>
298 return serializationTestObject_<General3rdOrderController>();
302template<
class TypeTag>
307 this->suggested_next_timestep_ = x;
310template<
class TypeTag>
315 return this->suggested_next_timestep_;
318template<
class TypeTag>
323 return *this->time_step_control_;
327template<
class TypeTag>
334 if (max_next_tstep > 0) {
335 this->suggested_next_timestep_ = max_next_tstep;
339template<
class TypeTag>
342updateTUNING(
double max_next_tstep,
const Tuning& tuning)
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;
348 updateNEXTSTEP(max_next_tstep);
349 this->timestep_after_event_ = tuning.TMAXWC;
356template<
class TypeTag>
357template<
class T,
class Serializer>
362 if (!serializer.isSerializing()) {
363 this->time_step_control_ = std::make_unique<T>();
365 serializer(*
static_cast<T*
>(this->time_step_control_.get()));
368template<
class TypeTag>
371AdaptiveTimeStepping<TypeTag>::
372castAndComp(
const AdaptiveTimeStepping<TypeTag>& Rhs)
const
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>::
382maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
const double original_time_step,
386 if (this->suggested_next_timestep_ < 0) {
387 this->suggested_next_timestep_ = this->restart_factor_ * original_time_step;
390 if (this->full_timestep_initially_) {
391 this->suggested_next_timestep_ = original_time_step;
395 if (is_event && this->timestep_after_event_ > 0) {
396 this->suggested_next_timestep_ = this->timestep_after_event_;
400template<
class TypeTag>
401template<
class Controller>
402AdaptiveTimeStepping<TypeTag>
403AdaptiveTimeStepping<TypeTag>::
404serializationTestObject_()
406 AdaptiveTimeStepping<TypeTag> result;
408 result.restart_factor_ = 1.0;
409 result.growth_factor_ = 2.0;
410 result.max_growth_ = 3.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>
434init_(
const UnitSystem& unitSystem)
436 std::tie(time_step_control_type_,
440 if (this->growth_factor_ < 1.0) {
441 OPM_THROW(std::runtime_error,
442 "Growth factor cannot be less than 1.");
452template<
class TypeTag>
453template<
class Solver>
459 const TuningUpdateCallback& tuning_updater)
460 : adaptive_time_stepping_{adaptive_time_stepping}
461 , simulator_timer_{simulator_timer}
463 , is_event_{is_event}
464 , tuning_updater_{tuning_updater}
468template<
class TypeTag>
469template<
class Solver>
470AdaptiveTimeStepping<TypeTag>&
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
484 if (isReservoirCouplingSlave_() && reservoirCouplingSlave_().activated()) {
485 return runStepReservoirCouplingSlave_();
487 else if (isReservoirCouplingMaster_() && reservoirCouplingMaster_().activated()) {
488 return runStepReservoirCouplingMaster_();
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>::
526maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
const double original_time_step)
528 this->adaptive_time_stepping_.maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
529 original_time_step, this->is_event_
536template<
class TypeTag>
537template<
class Solver>
539AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
540maybeUpdateTuning_(
double elapsed,
double dt,
int sub_step_number)
const
542 return this->tuning_updater_(elapsed, dt, sub_step_number);
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();
561 const auto original_time_step = this->simulator_timer_.currentStepLength();
562 const auto report_step = this->simulator_timer_.reportStepNum();
563 maybeUpdateTuning_(elapsed, original_time_step, report_step);
564 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(original_time_step);
566 AdaptiveSimulatorTimer substep_timer{
567 this->simulator_timer_.startDateTime(),
570 suggestedNextTimestep_(),
574 SubStepIteration<Solver> substepIteration{*
this, substep_timer, original_time_step,
true};
575 return substepIteration.run();
578#ifdef RESERVOIR_COUPLING_ENABLED
579template <
class TypeTag>
580template <
class Solver>
581ReservoirCouplingMaster<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
582AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
583reservoirCouplingMaster_()
585 return *(this->solver_.model().simulator().reservoirCouplingMaster());
589#ifdef RESERVOIR_COUPLING_ENABLED
590template <
class TypeTag>
591template <
class Solver>
592ReservoirCouplingSlave<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
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_()
638 const double original_time_step = this->simulator_timer_.currentStepLength();
639 double current_time{this->simulator_timer_.simulationTimeElapsed()};
640 double step_end_time = current_time + original_time_step;
641 auto current_step_length = original_time_step;
642 auto report_step_idx = this->simulator_timer_.currentStepNum();
643 if (report_step_idx == 0 && iteration == 0) {
644 reservoirCouplingMaster_().initTimeStepping();
646 SimulatorReport report;
648 reservoirCouplingMaster_().maybeReceiveActivationHandshakeFromSlaves(current_time);
650 reservoirCouplingMaster_().receiveNextReportDateFromSlaves();
651 bool start_of_report_step = (iteration == 0);
652 if (start_of_report_step) {
653 reservoirCouplingMaster_().initStartOfReportStep(report_step_idx);
655 current_step_length = reservoirCouplingMaster_().maybeChopSubStep(
656 current_step_length, current_time);
657 reservoirCouplingMaster_().sendNextTimeStepToSlaves(current_step_length);
658 if (start_of_report_step) {
659 maybeUpdateTuning_(current_time, current_step_length, 0);
660 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
662 AdaptiveSimulatorTimer substep_timer{
663 this->simulator_timer_.startDateTime(),
666 suggestedNextTimestep_(),
667 this->simulator_timer_.reportStepNum(),
671 current_time + current_step_length, step_end_time
673 SubStepIteration<Solver> substepIteration{*
this, substep_timer, current_step_length, final_step};
674 const auto sub_steps_report = substepIteration.run();
675 report += sub_steps_report;
676 current_time += current_step_length;
686#ifdef RESERVOIR_COUPLING_ENABLED
687template <
class TypeTag>
688template <
class Solver>
690AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
691runStepReservoirCouplingSlave_()
694 const double original_time_step = this->simulator_timer_.currentStepLength();
695 double current_time{this->simulator_timer_.simulationTimeElapsed()};
696 double step_end_time = current_time + original_time_step;
697 SimulatorReport report;
698 auto report_step_idx = this->simulator_timer_.currentStepNum();
699 if (report_step_idx == 0 && iteration == 0) {
700 reservoirCouplingSlave_().initTimeStepping();
703 bool start_of_report_step = (iteration == 0);
704 reservoirCouplingSlave_().sendNextReportDateToMasterProcess();
705 const auto timestep = reservoirCouplingSlave_().receiveNextTimeStepFromMaster();
706 if (start_of_report_step) {
707 maybeUpdateTuning_(current_time, original_time_step, 0);
708 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
710 AdaptiveSimulatorTimer substep_timer{
711 this->simulator_timer_.startDateTime(),
714 suggestedNextTimestep_(),
715 this->simulator_timer_.reportStepNum(),
719 current_time + timestep, step_end_time
721 SubStepIteration<Solver> substepIteration{*
this, substep_timer, timestep, final_step};
722 const auto sub_steps_report = substepIteration.run();
723 report += sub_steps_report;
724 current_time += timestep;
735template <
class TypeTag>
736template <
class Solver>
738AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
739suggestedNextTimestep_()
const
741 return this->adaptive_time_stepping_.suggestedNextStep();
750template<
class TypeTag>
751template<
class Solver>
752AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
753SubStepIteration(SubStepper<Solver>& substepper,
754 AdaptiveSimulatorTimer& substep_timer,
755 const double original_time_step,
757 : substepper_{substepper}
758 , substep_timer_{substep_timer}
759 , original_time_step_{original_time_step}
760 , final_step_{final_step}
761 , adaptive_time_stepping_{substepper.getAdaptiveTimerStepper()}
765template <
class TypeTag>
766template <
class Solver>
768AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
771 auto& simulator = solver_().model().simulator();
772 auto& problem = simulator.problem();
775 SimulatorReport report;
778 while (!this->substep_timer_.done()) {
782 maybeUpdateTuningAndTimeStep_();
784 const double dt = this->substep_timer_.currentStepLength();
785 if (timeStepVerbose_()) {
789 const auto substep_report = runSubStep_();
792 problem.setSubStepReport(substep_report);
793 auto& full_report = adaptive_time_stepping_.report();
794 full_report += substep_report;
795 problem.setSimulationReport(full_report);
797 report += substep_report;
799 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
800 ++this->substep_timer_;
802 const int iterations = getNumIterations_(substep_report);
803 auto dt_estimate = timeStepControlComputeEstimate_(
804 dt, iterations, this->substep_timer_);
806 assert(dt_estimate > 0);
807 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
810 maybeReportSubStep_(substep_report);
811 if (this->final_step_ && this->substep_timer_.done()) {
816 report.success.output_write_time += writeOutput_();
820 setTimeStep_(dt_estimate);
822 report.success.converged = this->substep_timer_.done();
823 this->substep_timer_.setLastStepFailed(
false);
826 this->substep_timer_.setLastStepFailed(
true);
827 checkTimeStepMaxRestartLimit_(restarts);
829 double new_time_step = restartFactor_() * dt;
830 if (substep_report.time_step_rejected) {
831 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>();
832 const double safetyFactor = Parameters::Get<Parameters::TimeStepControlSafetyFactor>();
833 const double temp_time_step = std::sqrt(safetyFactor * tol / solver_().model().relativeChange()) * dt;
834 if (temp_time_step < dt) {
835 new_time_step = temp_time_step;
838 checkTimeStepMinLimit_(new_time_step);
839 bool wells_shut =
false;
840 if (new_time_step > minTimeStepBeforeClosingWells_()) {
841 chopTimeStep_(new_time_step);
843 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
852 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
854 updateSuggestedNextStep_();
864template<
class TypeTag>
865template<
class Solver>
867AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
868checkContinueOnUnconvergedSolution_(
double dt)
const
870 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
871 if (continue_on_uncoverged_solution && solverVerbose_()) {
873 const auto msg = fmt::format(
874 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
875 "which is the minimum threshold given by option --solver-min-time-step\n",
878 OpmLog::problem(msg);
880 return continue_on_uncoverged_solution;
883template<
class TypeTag>
884template<
class Solver>
886AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
887checkTimeStepMaxRestartLimit_(
const int restarts)
const
891 if (restarts >= solverRestartMax_()) {
892 const auto msg = fmt::format(
893 "Solver failed to converge after cutting timestep {} times.", restarts
895 if (solverVerbose_()) {
899 throw TimeSteppingBreakdown{msg};
903template<
class TypeTag>
904template<
class Solver>
906AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
907checkTimeStepMinLimit_(
const double new_time_step)
const
909 using Meas = UnitSystem::measure;
912 if (new_time_step < minTimeStep_()) {
913 auto msg = fmt::format(
"Solver failed to converge after cutting timestep to ");
914 if (Parameters::Get<Parameters::EnableTuning>()) {
915 const UnitSystem& unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
917 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
918 unit_system.from_si(Meas::time, minTimeStep_()),
919 unit_system.name(Meas::time)
924 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
925 minTimeStep_() / 86400.0
928 if (solverVerbose_()) {
932 throw TimeSteppingBreakdown{msg};
936template<
class TypeTag>
937template<
class Solver>
939AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
940chopTimeStep_(
const double new_time_step)
942 setTimeStep_(new_time_step);
943 if (solverVerbose_()) {
944 const auto msg = fmt::format(
"{}\nTimestep chopped to {} days\n",
945 this->cause_of_failure_,
946 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
947 OpmLog::problem(msg);
951template<
class TypeTag>
952template<
class Solver>
954AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
955chopTimeStepOrCloseFailingWells_(
const double new_time_step)
957 bool wells_shut =
false;
964 const bool requireRepeatedFailures =
965 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
966 const std::set<std::string> failing_wells =
969 if (failing_wells.empty()) {
971 chopTimeStep_(new_time_step);
974 std::vector<std::string> shut_wells;
975 for (
const auto& well : failing_wells) {
976 const bool was_shut =
977 solver_().model().wellModel().forceShutWellByName(well,
978 this->substep_timer_.simulationTimeElapsed(),
981 shut_wells.push_back(well);
985 if (shut_wells.empty()) {
986 for (
const auto& well : failing_wells) {
987 const bool was_shut =
988 solver_().model().wellModel().forceShutWellByName(well,
989 this->substep_timer_.simulationTimeElapsed(),
992 shut_wells.push_back(well);
997 if (shut_wells.empty()) {
998 chopTimeStep_(new_time_step);
1001 if (solverVerbose_()) {
1002 const std::string msg =
1003 fmt::format(
"\nProblematic well(s) were shut: {}"
1004 "(retrying timestep)\n",
1005 fmt::join(shut_wells,
" "));
1006 OpmLog::problem(msg);
1013template<
class TypeTag>
1014template<
class Solver>
1015boost::posix_time::ptime
1016AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1017currentDateTime_()
const
1019 return simulatorTimer_().currentDateTime();
1022template<
class TypeTag>
1023template<
class Solver>
1025AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1026getNumIterations_(
const SimulatorReportSingle &substep_report)
const
1028 if (useNewtonIteration_()) {
1029 return substep_report.total_newton_iterations;
1032 return substep_report.total_linear_iterations;
1036template<
class TypeTag>
1037template<
class Solver>
1039AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1040growthFactor_()
const
1042 return this->adaptive_time_stepping_.growth_factor_;
1045template<
class TypeTag>
1046template<
class Solver>
1048AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1049ignoreConvergenceFailure_()
const
1051 return adaptive_time_stepping_.ignore_convergence_failure_;
1054template<
class TypeTag>
1055template<
class Solver>
1057AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1060 return this->adaptive_time_stepping_.max_growth_;
1063template<
class TypeTag>
1064template<
class Solver>
1066AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1067maybeReportSubStep_(SimulatorReportSingle substep_report)
const
1069 if (timeStepVerbose_()) {
1070 std::ostringstream ss;
1071 substep_report.reportStep(ss);
1072 OpmLog::info(ss.str());
1076template<
class TypeTag>
1077template<
class Solver>
1079AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1080maybeRestrictTimeStepGrowth_(
const double dt,
double dt_estimate,
const int restarts)
const
1083 dt_estimate = std::min(dt_estimate,
double(maxGrowth_() * dt));
1084 assert(dt_estimate > 0);
1087 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1096template<
class TypeTag>
1097template<
class Solver>
1099AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1100maybeUpdateTuningAndTimeStep_()
1108 const auto old_value = suggestedNextTimestep_();
1109 if (this->substepper_.maybeUpdateTuning_(this->substep_timer_.simulationTimeElapsed(),
1110 this->substep_timer_.currentStepLength(),
1111 this->substep_timer_.currentStepNum()))
1118 setTimeStep_(suggestedNextTimestep_());
1119 setSuggestedNextStep_(old_value);
1123template<
class TypeTag>
1124template<
class Solver>
1126AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1127minTimeStepBeforeClosingWells_()
const
1129 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1132template<
class TypeTag>
1133template<
class Solver>
1135AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1138 return this->adaptive_time_stepping_.min_time_step_;
1141template<
class TypeTag>
1142template<
class Solver>
1144AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1145restartFactor_()
const
1147 return this->adaptive_time_stepping_.restart_factor_;
1150template<
class TypeTag>
1151template<
class Solver>
1152SimulatorReportSingle
1153AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1156 SimulatorReportSingle substep_report;
1158 auto handleFailure = [
this, &substep_report]
1159 (
const std::string& failure_reason,
const std::exception& e,
bool log_exception =
true)
1161 substep_report = solver_().failureReport();
1162 this->cause_of_failure_ = failure_reason;
1163 if (log_exception && solverVerbose_()) {
1164 OpmLog::debug(std::string(
"Caught Exception: ") + e.what());
1169 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1170 if (solverVerbose_()) {
1172 OpmLog::debug(
"Overall linear iterations used: "
1176 catch (
const TooManyIterations& e) {
1177 handleFailure(
"Solver convergence failure - Iteration limit reached", e);
1179 catch (
const TimeSteppingBreakdown& e) {
1180 handleFailure(e.what(), e);
1182 catch (
const ConvergenceMonitorFailure& e) {
1183 handleFailure(
"Convergence monitor failure", e,
false);
1185 catch (
const LinearSolverProblem& e) {
1186 handleFailure(
"Linear solver convergence failure", e);
1188 catch (
const NumericalProblem& e) {
1189 handleFailure(
"Solver convergence failure - Numerical problem encountered", e);
1191 catch (
const std::runtime_error& e) {
1192 handleFailure(
"Runtime error encountered", e);
1194 catch (
const Dune::ISTLError& e) {
1195 handleFailure(
"ISTL error - Time step too large", e);
1197 catch (
const Dune::MatrixBlockError& e) {
1198 handleFailure(
"Matrix block error", e);
1201 return substep_report;
1204template<
class TypeTag>
1205template<
class Solver>
1207AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1208setTimeStep_(
double dt_estimate)
1210 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1213template<
class TypeTag>
1214template<
class Solver>
1216AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1219 return this->substepper_.solver_;
1223template<
class TypeTag>
1224template<
class Solver>
1226AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1227solverRestartMax_()
const
1229 return this->adaptive_time_stepping_.solver_restart_max_;
1232template<
class TypeTag>
1233template<
class Solver>
1235AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1236setSuggestedNextStep_(
double step)
1238 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1241template <
class TypeTag>
1242template <
class Solver>
1243const SimulatorTimer&
1244AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1245simulatorTimer_()
const
1247 return this->substepper_.simulator_timer_;
1250template <
class TypeTag>
1251template <
class Solver>
1253AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1254solverVerbose_()
const
1256 return this->adaptive_time_stepping_.solver_verbose_;
1259template<
class TypeTag>
1260template<
class Solver>
1261boost::posix_time::ptime
1262AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1263startDateTime_()
const
1265 return simulatorTimer_().startDateTime();
1268template <
class TypeTag>
1269template <
class Solver>
1271AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1272suggestedNextTimestep_()
const
1274 return this->adaptive_time_stepping_.suggestedNextStep();
1277template <
class TypeTag>
1278template <
class Solver>
1280AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1281timeStepControlComputeEstimate_(
const double dt,
const int iterations,
1282 const AdaptiveSimulatorTimer& substepTimer)
const
1285 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1286 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1287 dt, iterations, relative_change, substepTimer);
1290template <
class TypeTag>
1291template <
class Solver>
1293AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1294timeStepVerbose_()
const
1296 return this->adaptive_time_stepping_.timestep_verbose_;
1306template <
class TypeTag>
1307template <
class Solver>
1309AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1310updateSuggestedNextStep_()
1312 auto suggested_next_step = this->substep_timer_.currentStepLength();
1313 if (! std::isfinite(suggested_next_step)) {
1314 suggested_next_step = this->original_time_step_;
1316 if (timeStepVerbose_()) {
1317 std::ostringstream ss;
1318 this->substep_timer_.report(ss);
1319 ss <<
"Suggested next step size = "
1320 << unit::convert::to(suggested_next_step, unit::day) <<
" (days)" << std::endl;
1321 OpmLog::debug(ss.str());
1323 setSuggestedNextStep_(suggested_next_step);
1326template <
class TypeTag>
1327template <
class Solver>
1329AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1330useNewtonIteration_()
const
1332 return this->adaptive_time_stepping_.use_newton_iteration_;
1335template <
class TypeTag>
1336template <
class Solver>
1338AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1341 time::StopWatch perf_timer;
1343 auto& problem = solver_().model().simulator().problem();
1344 problem.writeOutput(
true);
1345 return perf_timer.secsSinceStart();
1352template<
class TypeTag>
1353template<
class Solver>
1354AdaptiveTimeStepping<TypeTag>::
1355SolutionTimeErrorSolverWrapper<Solver>::
1356SolutionTimeErrorSolverWrapper(
const Solver& solver)
1360template<
class TypeTag>
1361template<
class Solver>
1362double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange()
const
1365 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:434
void setSuggestedNextStep(const double x)
Definition: AdaptiveTimeStepping_impl.hpp:305
double suggestedNextStep() const
Definition: AdaptiveTimeStepping_impl.hpp:313
std::function< bool(const double, const double, const int)> TuningUpdateCallback
Definition: AdaptiveTimeStepping.hpp:80
bool operator==(const AdaptiveTimeStepping< TypeTag > &rhs) const
Definition: AdaptiveTimeStepping_impl.hpp:138
static AdaptiveTimeStepping< TypeTag > serializationTestObjectSimple()
Definition: AdaptiveTimeStepping_impl.hpp:288
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:217
void updateTUNING(double max_next_tstep, const Tuning &tuning)
Definition: AdaptiveTimeStepping_impl.hpp:342
TimeStepControlType time_step_control_type_
type of time step control object
Definition: AdaptiveTimeStepping.hpp:248
const TimeStepControlInterface & timeStepControl() const
Definition: AdaptiveTimeStepping_impl.hpp:321
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
void updateNEXTSTEP(double max_next_tstep)
Definition: AdaptiveTimeStepping_impl.hpp:330
static AdaptiveTimeStepping< TypeTag > serializationTestObjectHardcoded()
Definition: AdaptiveTimeStepping_impl.hpp:264
AdaptiveTimeStepping()=default
TimeStepController time_step_control_
time step control object
Definition: AdaptiveTimeStepping.hpp:249
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPIDIt()
Definition: AdaptiveTimeStepping_impl.hpp:280
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:296
SimulatorReport & report()
Definition: AdaptiveTimeStepping_impl.hpp:256
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPID()
Definition: AdaptiveTimeStepping_impl.hpp:272
static void registerParameters()
Definition: AdaptiveTimeStepping_impl.hpp:184
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:187
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: blackoilbioeffectsmodules.hh:43
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