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/common/timer.hh>
31#include <dune/istl/istlexception.hh>
33#include <opm/common/Exceptions.hpp>
34#include <opm/common/ErrorMacros.hpp>
35#include <opm/common/OpmLog/OpmLog.hpp>
37#include <opm/grid/utility/StopWatch.hpp>
39#include <opm/input/eclipse/Schedule/Tuning.hpp>
41#include <opm/input/eclipse/Units/Units.hpp>
42#include <opm/input/eclipse/Units/UnitSystem.hpp>
54#include <boost/date_time/posix_time/posix_time.hpp>
55#include <fmt/format.h>
56#include <fmt/ranges.h>
64template<
class TypeTag>
68 const double max_next_tstep,
69 const bool terminal_output
71 : time_step_control_{}
72 , restart_factor_{Parameters::
Get<Parameters::SolverRestartFactor<Scalar>>()}
73 , growth_factor_{Parameters::
Get<Parameters::SolverGrowthFactor<Scalar>>()}
74 , max_growth_{Parameters::
Get<Parameters::SolverMaxGrowth<Scalar>>()}
76 Parameters::
Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60}
78 unit_system.to_si(UnitSystem::measure::time,
79 Parameters::
Get<Parameters::SolverMinTimeStep<Scalar>>())}
80 , ignore_convergence_failure_{
81 Parameters::
Get<Parameters::SolverContinueOnConvergenceFailure>()}
82 , solver_restart_max_{Parameters::
Get<Parameters::SolverMaxRestarts>()}
83 , solver_verbose_{Parameters::
Get<Parameters::SolverVerbosity>() > 0 && terminal_output}
84 , timestep_verbose_{Parameters::
Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output}
85 , suggested_next_timestep_{
86 (max_next_tstep <= 0 ? Parameters::
Get<Parameters::InitialTimeStepInDays>()
87 : max_next_tstep) * 24 * 60 * 60}
88 , full_timestep_initially_{Parameters::
Get<Parameters::FullTimeStepInitially>()}
89 , timestep_after_event_{
90 Parameters::
Get<Parameters::TimeStepAfterEventInDays<Scalar>>() * 24 * 60 * 60}
91 , use_newton_iteration_{false}
92 , min_time_step_before_shutting_problematic_wells_{
93 Parameters::
Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
105template<
class TypeTag>
108 const Tuning& tuning,
109 const UnitSystem& unit_system,
111 const bool terminal_output
113 : time_step_control_{}
114 , restart_factor_{tuning.TSFCNV}
115 , growth_factor_{tuning.TFDIFF}
116 , max_growth_{tuning.TSFMAX}
117 , max_time_step_{tuning.TSMAXZ}
118 , min_time_step_{tuning.TSMINZ}
119 , ignore_convergence_failure_{true}
120 , solver_restart_max_{Parameters::
Get<Parameters::SolverMaxRestarts>()}
121 , solver_verbose_{Parameters::
Get<Parameters::SolverVerbosity>() > 0 && terminal_output}
122 , timestep_verbose_{Parameters::
Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output}
123 , suggested_next_timestep_{
124 max_next_tstep <= 0 ? Parameters::
Get<Parameters::InitialTimeStepInDays>() * 24 * 60 * 60
126 , full_timestep_initially_{Parameters::
Get<Parameters::FullTimeStepInitially>()}
127 , timestep_after_event_{tuning.TMAXWC}
128 , use_newton_iteration_{false}
129 , min_time_step_before_shutting_problematic_wells_{
130 Parameters::
Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
136template<
class TypeTag>
148 switch (this->time_step_control_type_) {
150 result = castAndComp<HardcodedTimeStepControl>(rhs);
153 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
156 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
159 result = castAndComp<PIDTimeStepControl>(rhs);
162 result = castAndComp<General3rdOrderController>(rhs);
178 this->min_time_step_before_shutting_problematic_wells_ ==
182template<
class TypeTag>
187 registerEclTimeSteppingParameters<Scalar>();
199template<
class TypeTag>
200template <
class Solver>
208 SubStepper<Solver> sub_stepper{
209 *
this, simulator_timer, solver, is_event, tuning_updater,
211 return sub_stepper.run();
214template<
class TypeTag>
215template<
class Serializer>
220 serializer(this->time_step_control_type_);
221 switch (this->time_step_control_type_) {
223 allocAndSerialize<HardcodedTimeStepControl>(serializer);
226 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
229 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
232 allocAndSerialize<PIDTimeStepControl>(serializer);
235 allocAndSerialize<General3rdOrderController>(serializer);
238 serializer(this->restart_factor_);
239 serializer(this->growth_factor_);
240 serializer(this->max_growth_);
241 serializer(this->max_time_step_);
242 serializer(this->min_time_step_);
243 serializer(this->ignore_convergence_failure_);
244 serializer(this->solver_restart_max_);
245 serializer(this->solver_verbose_);
246 serializer(this->timestep_verbose_);
247 serializer(this->suggested_next_timestep_);
248 serializer(this->full_timestep_initially_);
249 serializer(this->timestep_after_event_);
250 serializer(this->use_newton_iteration_);
251 serializer(this->min_time_step_before_shutting_problematic_wells_);
254template<
class TypeTag>
262template<
class TypeTag>
267 return serializationTestObject_<HardcodedTimeStepControl>();
270template<
class TypeTag>
275 return serializationTestObject_<PIDTimeStepControl>();
278template<
class TypeTag>
283 return serializationTestObject_<PIDAndIterationCountTimeStepControl>();
286template<
class TypeTag>
291 return serializationTestObject_<SimpleIterationCountTimeStepControl>();
294template<
class TypeTag>
299 return serializationTestObject_<General3rdOrderController>();
303template<
class TypeTag>
308 this->suggested_next_timestep_ = x;
311template<
class TypeTag>
316 return this->suggested_next_timestep_;
319template<
class TypeTag>
324 return *this->time_step_control_;
328template<
class TypeTag>
335 if (max_next_tstep > 0) {
336 this->suggested_next_timestep_ = max_next_tstep;
340template<
class TypeTag>
343updateTUNING(
double max_next_tstep,
const Tuning& tuning)
345 this->restart_factor_ = tuning.TSFCNV;
346 this->growth_factor_ = tuning.TFDIFF;
347 this->max_growth_ = tuning.TSFMAX;
348 this->max_time_step_ = tuning.TSMAXZ;
349 updateNEXTSTEP(max_next_tstep);
350 this->timestep_after_event_ = tuning.TMAXWC;
357template<
class TypeTag>
358template<
class T,
class Serializer>
363 if (!serializer.isSerializing()) {
364 this->time_step_control_ = std::make_unique<T>();
366 serializer(*
static_cast<T*
>(this->time_step_control_.get()));
369template<
class TypeTag>
372AdaptiveTimeStepping<TypeTag>::
373castAndComp(
const AdaptiveTimeStepping<TypeTag>& Rhs)
const
375 const T* lhs =
static_cast<const T*
>(this->time_step_control_.get());
376 const T* rhs =
static_cast<const T*
>(Rhs.time_step_control_.get());
380template<
class TypeTag>
382AdaptiveTimeStepping<TypeTag>::
383maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
const double original_time_step,
387 if (this->suggested_next_timestep_ < 0) {
388 this->suggested_next_timestep_ = this->restart_factor_ * original_time_step;
391 if (this->full_timestep_initially_) {
392 this->suggested_next_timestep_ = original_time_step;
396 if (is_event && this->timestep_after_event_ > 0) {
397 this->suggested_next_timestep_ = this->timestep_after_event_;
401template<
class TypeTag>
402template<
class Controller>
403AdaptiveTimeStepping<TypeTag>
404AdaptiveTimeStepping<TypeTag>::
405serializationTestObject_()
407 AdaptiveTimeStepping<TypeTag> result;
409 result.restart_factor_ = 1.0;
410 result.growth_factor_ = 2.0;
411 result.max_growth_ = 3.0;
412 result.max_time_step_ = 4.0;
413 result.min_time_step_ = 5.0;
414 result.ignore_convergence_failure_ =
true;
415 result.solver_restart_max_ = 6;
416 result.solver_verbose_ =
true;
417 result.timestep_verbose_ =
true;
418 result.suggested_next_timestep_ = 7.0;
419 result.full_timestep_initially_ =
true;
420 result.use_newton_iteration_ =
true;
421 result.min_time_step_before_shutting_problematic_wells_ = 9.0;
422 result.time_step_control_type_ = Controller::Type;
423 result.time_step_control_ =
424 std::make_unique<Controller>(Controller::serializationTestObject());
433template<
class TypeTag>
435init_(
const UnitSystem& unitSystem)
437 std::tie(time_step_control_type_,
441 if (this->growth_factor_ < 1.0) {
442 OPM_THROW(std::runtime_error,
443 "Growth factor cannot be less than 1.");
453template<
class TypeTag>
454template<
class Solver>
460 const TuningUpdateCallback& tuning_updater)
461 : adaptive_time_stepping_{adaptive_time_stepping}
462 , simulator_timer_{simulator_timer}
464 , is_event_{is_event}
465 , tuning_updater_{tuning_updater}
469template<
class TypeTag>
470template<
class Solver>
471AdaptiveTimeStepping<TypeTag>&
472AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
473getAdaptiveTimerStepper()
475 return adaptive_time_stepping_;
478template<
class TypeTag>
479template<
class Solver>
481AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
484#ifdef RESERVOIR_COUPLING_ENABLED
485 if (isReservoirCouplingSlave_() && reservoirCouplingSlave_().activated()) {
486 return runStepReservoirCouplingSlave_();
488 else if (isReservoirCouplingMaster_() && reservoirCouplingMaster_().activated()) {
489 return runStepReservoirCouplingMaster_();
492 return runStepOriginal_();
495 return runStepOriginal_();
503#ifdef RESERVOIR_COUPLING_ENABLED
504template<
class TypeTag>
505template<
class Solver>
507AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
508isReservoirCouplingMaster_()
const
510 return this->solver_.model().simulator().reservoirCouplingMaster() !=
nullptr;
513template<
class TypeTag>
514template<
class Solver>
516AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
517isReservoirCouplingSlave_()
const
519 return this->solver_.model().simulator().reservoirCouplingSlave() !=
nullptr;
523template<
class TypeTag>
524template<
class Solver>
526AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
527maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
const double original_time_step)
529 this->adaptive_time_stepping_.maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
530 original_time_step, this->is_event_
537template<
class TypeTag>
538template<
class Solver>
540AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
541maybeUpdateTuning_(
double elapsed,
double dt,
int sub_step_number)
const
543 return this->tuning_updater_(elapsed, dt, sub_step_number);
546template<
class TypeTag>
547template<
class Solver>
549AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
552 return this->adaptive_time_stepping_.max_time_step_;
555template <
class TypeTag>
556template <
class Solver>
558AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
561 const auto elapsed = this->simulator_timer_.simulationTimeElapsed();
562 const auto original_time_step = this->simulator_timer_.currentStepLength();
563 const auto report_step = this->simulator_timer_.reportStepNum();
564 maybeUpdateTuning_(elapsed, original_time_step, report_step);
565 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(original_time_step);
567 AdaptiveSimulatorTimer substep_timer{
568 this->simulator_timer_.startDateTime(),
571 suggestedNextTimestep_(),
575 SubStepIteration<Solver> substepIteration{*
this, substep_timer, original_time_step,
true};
576 return substepIteration.run();
579#ifdef RESERVOIR_COUPLING_ENABLED
580template <
class TypeTag>
581template <
class Solver>
582ReservoirCouplingMaster<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
583AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
584reservoirCouplingMaster_()
586 return *(this->solver_.model().simulator().reservoirCouplingMaster());
590#ifdef RESERVOIR_COUPLING_ENABLED
591template <
class TypeTag>
592template <
class Solver>
593ReservoirCouplingSlave<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
594AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
595reservoirCouplingSlave_()
597 return *(this->solver_.model().simulator().reservoirCouplingSlave());
601#ifdef RESERVOIR_COUPLING_ENABLED
632template <
class TypeTag>
633template <
class Solver>
635AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
636runStepReservoirCouplingMaster_()
639 const double original_time_step = this->simulator_timer_.currentStepLength();
640 double current_time{this->simulator_timer_.simulationTimeElapsed()};
641 double step_end_time = current_time + original_time_step;
642 auto current_step_length = original_time_step;
643 auto report_step_idx = this->simulator_timer_.currentStepNum();
644 if (report_step_idx == 0 && iteration == 0) {
645 reservoirCouplingMaster_().initTimeStepping();
647 SimulatorReport report;
649 reservoirCouplingMaster_().maybeReceiveActivationHandshakeFromSlaves(current_time);
651 reservoirCouplingMaster_().sendDontTerminateSignalToSlaves();
652 reservoirCouplingMaster_().receiveNextReportDateFromSlaves();
653 bool start_of_report_step = (iteration == 0);
654 if (start_of_report_step) {
655 reservoirCouplingMaster_().initStartOfReportStep(report_step_idx);
657 current_step_length = reservoirCouplingMaster_().maybeChopSubStep(
658 current_step_length, current_time);
659 auto num_active = reservoirCouplingMaster_().numActivatedSlaves();
660 OpmLog::info(fmt::format(
661 "\nChoosing next sync time between master and {} active slave {}: {:.2f} days",
662 num_active, (num_active == 1 ?
"process" :
"processes"),
663 current_step_length / unit::day
665 reservoirCouplingMaster_().sendNextTimeStepToSlaves(current_step_length);
666 if (start_of_report_step) {
667 maybeUpdateTuning_(current_time, current_step_length, 0);
668 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
670 AdaptiveSimulatorTimer substep_timer{
671 this->simulator_timer_.startDateTime(),
674 suggestedNextTimestep_(),
675 this->simulator_timer_.reportStepNum(),
679 current_time + current_step_length, step_end_time
684 reservoirCouplingMaster_().setFirstSubstepOfSyncTimestep(
true);
685 SubStepIteration<Solver> substepIteration{*
this, substep_timer, current_step_length, final_step};
686 const auto sub_steps_report = substepIteration.run();
687 report += sub_steps_report;
688 current_time += current_step_length;
698#ifdef RESERVOIR_COUPLING_ENABLED
699template <
class TypeTag>
700template <
class Solver>
702AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
703runStepReservoirCouplingSlave_()
706 const double original_time_step = this->simulator_timer_.currentStepLength();
707 double current_time{this->simulator_timer_.simulationTimeElapsed()};
708 double step_end_time = current_time + original_time_step;
709 SimulatorReport report;
710 auto report_step_idx = this->simulator_timer_.currentStepNum();
711 if (report_step_idx == 0 && iteration == 0) {
712 reservoirCouplingSlave_().initTimeStepping();
715 bool start_of_report_step = (iteration == 0);
716 if (reservoirCouplingSlave_().maybeReceiveTerminateSignalFromMaster()) {
720 reservoirCouplingSlave_().sendNextReportDateToMasterProcess();
721 const auto timestep = reservoirCouplingSlave_().receiveNextTimeStepFromMaster();
722 if (start_of_report_step) {
723 maybeUpdateTuning_(current_time, original_time_step, 0);
724 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
726 AdaptiveSimulatorTimer substep_timer{
727 this->simulator_timer_.startDateTime(),
730 suggestedNextTimestep_(),
731 this->simulator_timer_.reportStepNum(),
735 current_time + timestep, step_end_time
740 reservoirCouplingSlave_().setFirstSubstepOfSyncTimestep(
true);
741 SubStepIteration<Solver> substepIteration{*
this, substep_timer, timestep, final_step};
742 const auto sub_steps_report = substepIteration.run();
743 report += sub_steps_report;
744 current_time += timestep;
755template <
class TypeTag>
756template <
class Solver>
758AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
759suggestedNextTimestep_()
const
761 return this->adaptive_time_stepping_.suggestedNextStep();
770template<
class TypeTag>
771template<
class Solver>
772AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
773SubStepIteration(SubStepper<Solver>& substepper,
774 AdaptiveSimulatorTimer& substep_timer,
775 const double original_time_step,
777 : substepper_{substepper}
778 , substep_timer_{substep_timer}
779 , original_time_step_{original_time_step}
780 , final_step_{final_step}
781 , adaptive_time_stepping_{substepper.getAdaptiveTimerStepper()}
785template <
class TypeTag>
786template <
class Solver>
788AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
791 auto& simulator = solver_().model().simulator();
792 auto& problem = simulator.problem();
795 SimulatorReport report;
798 while (!this->substep_timer_.done()) {
802 maybeUpdateTuningAndTimeStep_();
804 const double dt = this->substep_timer_.currentStepLength();
805 if (timeStepVerbose_()) {
809 auto substep_report = runSubStep_();
810 markFirstSubStepAsFinished_();
812 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
813 Dune::Timer perfTimer;
816 problem.setSubStepReport(substep_report);
817 auto& full_report = adaptive_time_stepping_.report();
818 full_report += substep_report;
819 problem.setSimulationReport(full_report);
820 problem.endTimeStep();
821 substep_report.pre_post_time += perfTimer.stop();
823 report += substep_report;
825 ++this->substep_timer_;
827 const int iterations = getNumIterations_(substep_report);
828 auto dt_estimate = timeStepControlComputeEstimate_(
829 dt, iterations, this->substep_timer_);
831 assert(dt_estimate > 0);
832 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
835 maybeReportSubStep_(substep_report);
836 if (this->final_step_ && this->substep_timer_.done()) {
841 report.success.output_write_time += writeOutput_();
845 setTimeStep_(dt_estimate);
847 report.success.converged = this->substep_timer_.done();
848 this->substep_timer_.setLastStepFailed(
false);
851 report += substep_report;
852 this->substep_timer_.setLastStepFailed(
true);
853 checkTimeStepMaxRestartLimit_(restarts);
855 double new_time_step = restartFactor_() * dt;
856 if (substep_report.time_step_rejected) {
857 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>();
858 const double safetyFactor = Parameters::Get<Parameters::TimeStepControlSafetyFactor>();
859 const double temp_time_step = std::sqrt(safetyFactor * tol / solver_().model().relativeChange()) * dt;
860 if (temp_time_step < dt) {
861 new_time_step = temp_time_step;
864 checkTimeStepMinLimit_(new_time_step);
865 bool wells_shut =
false;
866 if (new_time_step > minTimeStepBeforeClosingWells_()) {
867 chopTimeStep_(new_time_step);
869 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
878 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
880 updateSuggestedNextStep_();
890template<
class TypeTag>
891template<
class Solver>
893AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
894checkContinueOnUnconvergedSolution_(
double dt)
const
896 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
897 if (continue_on_uncoverged_solution && solverVerbose_()) {
899 const auto msg = fmt::format(
900 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
901 "which is the minimum threshold given by option --solver-min-time-step\n",
904 OpmLog::problem(msg);
906 return continue_on_uncoverged_solution;
909template<
class TypeTag>
910template<
class Solver>
912AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
913checkTimeStepMaxRestartLimit_(
const int restarts)
const
917 if (restarts >= solverRestartMax_()) {
918 const auto msg = fmt::format(
919 fmt::runtime(
"Solver failed to converge after cutting timestep {} times."), restarts
921 if (solverVerbose_()) {
925 throw TimeSteppingBreakdown{msg};
929template<
class TypeTag>
930template<
class Solver>
932AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
933checkTimeStepMinLimit_(
const double new_time_step)
const
935 using Meas = UnitSystem::measure;
938 if (new_time_step < minTimeStep_()) {
939 std::string msg =
"Solver failed to converge after cutting timestep to ";
940 if (Parameters::Get<Parameters::EnableTuning>()) {
941 const UnitSystem& unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
943 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
944 unit_system.from_si(Meas::time, minTimeStep_()),
945 unit_system.name(Meas::time)
950 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
951 minTimeStep_() / 86400.0
954 if (solverVerbose_()) {
958 throw TimeSteppingBreakdown{msg};
962template<
class TypeTag>
963template<
class Solver>
965AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
966chopTimeStep_(
const double new_time_step)
968 setTimeStep_(new_time_step);
969 if (solverVerbose_()) {
970 const auto msg = fmt::format(fmt::runtime(
"{}\nTimestep chopped to {} days\n"),
971 this->cause_of_failure_,
972 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
973 OpmLog::problem(msg);
977template<
class TypeTag>
978template<
class Solver>
980AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
981chopTimeStepOrCloseFailingWells_(
const double new_time_step)
983 bool wells_shut =
false;
990 const bool requireRepeatedFailures =
991 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
992 const std::set<std::string> failing_wells =
995 if (failing_wells.empty()) {
997 chopTimeStep_(new_time_step);
1000 std::vector<std::string> shut_wells;
1001 for (
const auto& well : failing_wells) {
1002 const bool was_shut =
1003 solver_().model().wellModel().forceShutWellByName(well,
1004 this->substep_timer_.simulationTimeElapsed(),
1007 shut_wells.push_back(well);
1011 if (shut_wells.empty()) {
1012 for (
const auto& well : failing_wells) {
1013 const bool was_shut =
1014 solver_().model().wellModel().forceShutWellByName(well,
1015 this->substep_timer_.simulationTimeElapsed(),
1018 shut_wells.push_back(well);
1023 if (shut_wells.empty()) {
1024 chopTimeStep_(new_time_step);
1027 if (solverVerbose_()) {
1028 const std::string msg =
1029 fmt::format(fmt::runtime(
"\nProblematic well(s) were shut: {}"
1030 "(retrying timestep)\n"),
1031 fmt::join(shut_wells,
" "));
1032 OpmLog::problem(msg);
1039template<
class TypeTag>
1040template<
class Solver>
1041boost::posix_time::ptime
1042AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1043currentDateTime_()
const
1045 return simulatorTimer_().currentDateTime();
1048template<
class TypeTag>
1049template<
class Solver>
1051AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1052getNumIterations_(
const SimulatorReportSingle &substep_report)
const
1054 if (useNewtonIteration_()) {
1055 return substep_report.total_newton_iterations;
1058 return substep_report.total_linear_iterations;
1062template<
class TypeTag>
1063template<
class Solver>
1065AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1066growthFactor_()
const
1068 return this->adaptive_time_stepping_.growth_factor_;
1071template<
class TypeTag>
1072template<
class Solver>
1074AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1075ignoreConvergenceFailure_()
const
1077 return adaptive_time_stepping_.ignore_convergence_failure_;
1080template<
class TypeTag>
1081template<
class Solver>
1083AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1084isReservoirCouplingMaster_()
const
1086 return this->substepper_.isReservoirCouplingMaster_();
1089template<
class TypeTag>
1090template<
class Solver>
1092AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1093isReservoirCouplingSlave_()
const
1095 return this->substepper_.isReservoirCouplingSlave_();
1098template<
class TypeTag>
1099template<
class Solver>
1101AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1102markFirstSubStepAsFinished_()
const
1104#ifdef RESERVOIR_COUPLING_ENABLED
1108 if (isReservoirCouplingMaster_()) {
1109 reservoirCouplingMaster_().setFirstSubstepOfSyncTimestep(
false);
1111 else if (isReservoirCouplingSlave_()) {
1112 reservoirCouplingSlave_().setFirstSubstepOfSyncTimestep(
false);
1118template<
class TypeTag>
1119template<
class Solver>
1121AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1124 return this->adaptive_time_stepping_.max_growth_;
1127template<
class TypeTag>
1128template<
class Solver>
1130AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1131maybeReportSubStep_(SimulatorReportSingle substep_report)
const
1133 if (timeStepVerbose_()) {
1134 std::ostringstream ss;
1135 substep_report.reportStep(ss);
1136 OpmLog::info(ss.str());
1140template<
class TypeTag>
1141template<
class Solver>
1143AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1144maybeRestrictTimeStepGrowth_(
const double dt,
double dt_estimate,
const int restarts)
const
1147 dt_estimate = std::min(dt_estimate,
double(maxGrowth_() * dt));
1148 assert(dt_estimate > 0);
1151 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1160template<
class TypeTag>
1161template<
class Solver>
1163AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1164maybeUpdateTuningAndTimeStep_()
1172 const auto old_value = suggestedNextTimestep_();
1173 if (this->substepper_.maybeUpdateTuning_(this->substep_timer_.simulationTimeElapsed(),
1174 this->substep_timer_.currentStepLength(),
1175 this->substep_timer_.currentStepNum()))
1182 setTimeStep_(suggestedNextTimestep_());
1183 setSuggestedNextStep_(old_value);
1187template<
class TypeTag>
1188template<
class Solver>
1190AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1191minTimeStepBeforeClosingWells_()
const
1193 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1196template<
class TypeTag>
1197template<
class Solver>
1199AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1202 return this->adaptive_time_stepping_.min_time_step_;
1205#ifdef RESERVOIR_COUPLING_ENABLED
1206template<
class TypeTag>
1207template<
class Solver>
1208ReservoirCouplingMaster<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
1209AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1210reservoirCouplingMaster_()
const
1212 return this->substepper_.reservoirCouplingMaster_();
1215template<
class TypeTag>
1216template<
class Solver>
1217ReservoirCouplingSlave<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
1218AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1219reservoirCouplingSlave_()
const
1221 return this->substepper_.reservoirCouplingSlave_();
1225template<
class TypeTag>
1226template<
class Solver>
1228AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1229restartFactor_()
const
1231 return this->adaptive_time_stepping_.restart_factor_;
1234template<
class TypeTag>
1235template<
class Solver>
1236SimulatorReportSingle
1237AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1240 SimulatorReportSingle substep_report;
1242 auto handleFailure = [
this, &substep_report]
1243 (
const std::string& failure_reason,
const std::exception& e,
bool log_exception =
true)
1245 substep_report = solver_().failureReport();
1246 this->cause_of_failure_ = failure_reason;
1247 if (log_exception && solverVerbose_()) {
1248 OpmLog::debug(std::string(
"Caught Exception: ") + e.what());
1253 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1254 if (solverVerbose_()) {
1256 OpmLog::debug(
"Overall linear iterations used: "
1260 catch (
const TooManyIterations& e) {
1261 handleFailure(
"Solver convergence failure - Iteration limit reached", e);
1263 catch (
const TimeSteppingBreakdown& e) {
1264 handleFailure(e.what(), e);
1266 catch (
const ConvergenceMonitorFailure& e) {
1267 handleFailure(
"Convergence monitor failure", e,
false);
1269 catch (
const LinearSolverProblem& e) {
1270 handleFailure(
"Linear solver convergence failure", e);
1272 catch (
const NumericalProblem& e) {
1273 handleFailure(
"Solver convergence failure - Numerical problem encountered", e);
1275 catch (
const std::runtime_error& e) {
1276 handleFailure(
"Runtime error encountered", e);
1278 catch (
const Dune::ISTLError& e) {
1279 handleFailure(
"ISTL error - Time step too large", e);
1281 catch (
const Dune::MatrixBlockError& e) {
1282 handleFailure(
"Matrix block error", e);
1285 return substep_report;
1288template<
class TypeTag>
1289template<
class Solver>
1291AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1292setTimeStep_(
double dt_estimate)
1294 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1297template<
class TypeTag>
1298template<
class Solver>
1300AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1303 return this->substepper_.solver_;
1307template<
class TypeTag>
1308template<
class Solver>
1310AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1311solverRestartMax_()
const
1313 return this->adaptive_time_stepping_.solver_restart_max_;
1316template<
class TypeTag>
1317template<
class Solver>
1319AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1320setSuggestedNextStep_(
double step)
1322 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1325template <
class TypeTag>
1326template <
class Solver>
1327const SimulatorTimer&
1328AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1329simulatorTimer_()
const
1331 return this->substepper_.simulator_timer_;
1334template <
class TypeTag>
1335template <
class Solver>
1337AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1338solverVerbose_()
const
1340 return this->adaptive_time_stepping_.solver_verbose_;
1343template<
class TypeTag>
1344template<
class Solver>
1345boost::posix_time::ptime
1346AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1347startDateTime_()
const
1349 return simulatorTimer_().startDateTime();
1352template <
class TypeTag>
1353template <
class Solver>
1355AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1356suggestedNextTimestep_()
const
1358 return this->adaptive_time_stepping_.suggestedNextStep();
1361template <
class TypeTag>
1362template <
class Solver>
1364AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1365timeStepControlComputeEstimate_(
const double dt,
const int iterations,
1366 const AdaptiveSimulatorTimer& substepTimer)
const
1369 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1370 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1371 dt, iterations, relative_change, substepTimer);
1374template <
class TypeTag>
1375template <
class Solver>
1377AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1378timeStepVerbose_()
const
1380 return this->adaptive_time_stepping_.timestep_verbose_;
1390template <
class TypeTag>
1391template <
class Solver>
1393AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1394updateSuggestedNextStep_()
1396 auto suggested_next_step = this->substep_timer_.currentStepLength();
1397 if (! std::isfinite(suggested_next_step)) {
1398 suggested_next_step = this->original_time_step_;
1400 if (timeStepVerbose_()) {
1401 std::ostringstream ss;
1402 this->substep_timer_.report(ss);
1403 ss <<
"Suggested next step size = "
1404 << unit::convert::to(suggested_next_step, unit::day) <<
" (days)" << std::endl;
1405 OpmLog::debug(ss.str());
1407 setSuggestedNextStep_(suggested_next_step);
1410template <
class TypeTag>
1411template <
class Solver>
1413AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1414useNewtonIteration_()
const
1416 return this->adaptive_time_stepping_.use_newton_iteration_;
1419template <
class TypeTag>
1420template <
class Solver>
1422AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1425 time::StopWatch perf_timer;
1427 auto& problem = solver_().model().simulator().problem();
1428 problem.writeOutput(
true);
1429 return perf_timer.secsSinceStart();
1436template<
class TypeTag>
1437template<
class Solver>
1438AdaptiveTimeStepping<TypeTag>::
1439SolutionTimeErrorSolverWrapper<Solver>::
1440SolutionTimeErrorSolverWrapper(
const Solver& solver)
1444template<
class TypeTag>
1445template<
class Solver>
1446double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange()
const
1449 return solver_.model().relativeChange();
Definition: AdaptiveTimeStepping.hpp:76
double max_growth_
factor that limits the maximum growth of a time step
Definition: AdaptiveTimeStepping.hpp:257
double max_time_step_
maximal allowed time step size in days
Definition: AdaptiveTimeStepping.hpp:258
bool solver_verbose_
solver verbosity
Definition: AdaptiveTimeStepping.hpp:262
int solver_restart_max_
how many restart of solver are allowed
Definition: AdaptiveTimeStepping.hpp:261
double timestep_after_event_
suggested size of timestep after an event
Definition: AdaptiveTimeStepping.hpp:266
void init_(const UnitSystem &unitSystem)
Definition: AdaptiveTimeStepping_impl.hpp:435
void setSuggestedNextStep(const double x)
Definition: AdaptiveTimeStepping_impl.hpp:306
double suggestedNextStep() const
Definition: AdaptiveTimeStepping_impl.hpp:314
std::function< bool(const double, const double, const int)> TuningUpdateCallback
Definition: AdaptiveTimeStepping.hpp:78
bool operator==(const AdaptiveTimeStepping< TypeTag > &rhs) const
Definition: AdaptiveTimeStepping_impl.hpp:139
static AdaptiveTimeStepping< TypeTag > serializationTestObjectSimple()
Definition: AdaptiveTimeStepping_impl.hpp:289
bool ignore_convergence_failure_
continue instead of stop when minimum time step is reached
Definition: AdaptiveTimeStepping.hpp:260
void serializeOp(Serializer &serializer)
Definition: AdaptiveTimeStepping_impl.hpp:218
void updateTUNING(double max_next_tstep, const Tuning &tuning)
Definition: AdaptiveTimeStepping_impl.hpp:343
TimeStepControlType time_step_control_type_
type of time step control object
Definition: AdaptiveTimeStepping.hpp:253
const TimeStepControlInterface & timeStepControl() const
Definition: AdaptiveTimeStepping_impl.hpp:322
bool full_timestep_initially_
beginning with the size of the time step from data file
Definition: AdaptiveTimeStepping.hpp:265
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:203
double growth_factor_
factor to multiply time step when solver recovered from failed convergence
Definition: AdaptiveTimeStepping.hpp:256
double restart_factor_
factor to multiply time step with when solver fails to converge
Definition: AdaptiveTimeStepping.hpp:255
double min_time_step_
minimal allowed time step size before throwing
Definition: AdaptiveTimeStepping.hpp:259
void updateNEXTSTEP(double max_next_tstep)
Definition: AdaptiveTimeStepping_impl.hpp:331
static AdaptiveTimeStepping< TypeTag > serializationTestObjectHardcoded()
Definition: AdaptiveTimeStepping_impl.hpp:265
AdaptiveTimeStepping()=default
TimeStepController time_step_control_
time step control object
Definition: AdaptiveTimeStepping.hpp:254
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPIDIt()
Definition: AdaptiveTimeStepping_impl.hpp:281
double min_time_step_before_shutting_problematic_wells_
< shut problematic wells when time step size in days are less than this
Definition: AdaptiveTimeStepping.hpp:270
static AdaptiveTimeStepping< TypeTag > serializationTestObject3rdOrder()
Definition: AdaptiveTimeStepping_impl.hpp:297
SimulatorReport & report()
Definition: AdaptiveTimeStepping_impl.hpp:257
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPID()
Definition: AdaptiveTimeStepping_impl.hpp:273
static void registerParameters()
Definition: AdaptiveTimeStepping_impl.hpp:185
bool use_newton_iteration_
use newton iteration count for adaptive time step control
Definition: AdaptiveTimeStepping.hpp:267
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