AdaptiveTimeStepping_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2024 Equinor ASA.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
21#define OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
22
23// Improve IDE experience
24#ifndef OPM_ADAPTIVE_TIME_STEPPING_HPP
25#include <config.h>
28#endif
29
30#include <dune/istl/istlexception.hh>
31
32#include <opm/common/Exceptions.hpp>
33#include <opm/common/ErrorMacros.hpp>
34#include <opm/common/OpmLog/OpmLog.hpp>
35
36#include <opm/grid/utility/StopWatch.hpp>
37
38#include <opm/input/eclipse/Schedule/Tuning.hpp>
39
40#include <opm/input/eclipse/Units/Units.hpp>
41#include <opm/input/eclipse/Units/UnitSystem.hpp>
42
44
46
47#include <algorithm>
48#include <cassert>
49#include <cmath>
50#include <sstream>
51#include <stdexcept>
52
53#include <boost/date_time/posix_time/posix_time.hpp>
54#include <fmt/format.h>
55#include <fmt/ranges.h>
56namespace Opm {
57/*********************************************
58 * Public methods of AdaptiveTimeStepping
59 * ******************************************/
60
61
63template<class TypeTag>
65AdaptiveTimeStepping(const UnitSystem& unit_system,
66 const SimulatorReport& report,
67 const double max_next_tstep,
68 const bool terminal_output
69)
70 : time_step_control_{}
71 , restart_factor_{Parameters::Get<Parameters::SolverRestartFactor<Scalar>>()} // 0.33
72 , growth_factor_{Parameters::Get<Parameters::SolverGrowthFactor<Scalar>>()} // 2.0
73 , max_growth_{Parameters::Get<Parameters::SolverMaxGrowth<Scalar>>()} // 3.0
74 , max_time_step_{
75 Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60} // 365.25
76 , min_time_step_{
77 unit_system.to_si(UnitSystem::measure::time,
78 Parameters::Get<Parameters::SolverMinTimeStep<Scalar>>())} // 1e-12;
79 , ignore_convergence_failure_{
80 Parameters::Get<Parameters::SolverContinueOnConvergenceFailure>()} // false;
81 , solver_restart_max_{Parameters::Get<Parameters::SolverMaxRestarts>()} // 10
82 , solver_verbose_{Parameters::Get<Parameters::SolverVerbosity>() > 0 && terminal_output} // 2
83 , timestep_verbose_{Parameters::Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output} // 2
84 , suggested_next_timestep_{
85 (max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>()
86 : max_next_tstep) * 24 * 60 * 60} // 1.0
87 , full_timestep_initially_{Parameters::Get<Parameters::FullTimeStepInitially>()} // false
88 , timestep_after_event_{
89 Parameters::Get<Parameters::TimeStepAfterEventInDays<Scalar>>() * 24 * 60 * 60} // 1e30
90 , use_newton_iteration_{false}
91 , min_time_step_before_shutting_problematic_wells_{
92 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
93 , report_(report)
94{
95 init_(unit_system);
96}
97
104template<class TypeTag>
106AdaptiveTimeStepping(double max_next_tstep,
107 const Tuning& tuning,
108 const UnitSystem& unit_system,
109 const SimulatorReport& report,
110 const bool terminal_output
111)
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} // 365.25
117 , min_time_step_{tuning.TSMINZ} // 0.1;
118 , ignore_convergence_failure_{true}
119 , solver_restart_max_{Parameters::Get<Parameters::SolverMaxRestarts>()} // 10
120 , solver_verbose_{Parameters::Get<Parameters::SolverVerbosity>() > 0 && terminal_output} // 2
121 , timestep_verbose_{Parameters::Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output} // 2
122 , suggested_next_timestep_{
123 max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>() * 24 * 60 * 60
124 : max_next_tstep} // 1.0
125 , full_timestep_initially_{Parameters::Get<Parameters::FullTimeStepInitially>()} // false
126 , timestep_after_event_{tuning.TMAXWC} // 1e30
127 , use_newton_iteration_{false}
128 , min_time_step_before_shutting_problematic_wells_{
129 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
130 , report_(report)
131{
132 init_(unit_system);
133}
134
135template<class TypeTag>
136bool
139{
140 if (this->time_step_control_type_ != rhs.time_step_control_type_ ||
141 (this->time_step_control_ && !rhs.time_step_control_) ||
142 (!this->time_step_control_ && rhs.time_step_control_)) {
143 return false;
144 }
145
146 bool result = false;
147 switch (this->time_step_control_type_) {
149 result = castAndComp<HardcodedTimeStepControl>(rhs);
150 break;
152 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
153 break;
155 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
156 break;
158 result = castAndComp<PIDTimeStepControl>(rhs);
159 break;
161 result = castAndComp<General3rdOrderController>(rhs);
162 break;
163 }
164
165 return result &&
166 this->restart_factor_ == rhs.restart_factor_ &&
167 this->growth_factor_ == rhs.growth_factor_ &&
168 this->max_growth_ == rhs.max_growth_ &&
169 this->max_time_step_ == rhs.max_time_step_ &&
170 this->min_time_step_ == rhs.min_time_step_ &&
171 this->ignore_convergence_failure_ == rhs.ignore_convergence_failure_ &&
172 this->solver_restart_max_== rhs.solver_restart_max_ &&
173 this->solver_verbose_ == rhs.solver_verbose_ &&
174 this->full_timestep_initially_ == rhs.full_timestep_initially_ &&
175 this->timestep_after_event_ == rhs.timestep_after_event_ &&
176 this->use_newton_iteration_ == rhs.use_newton_iteration_ &&
177 this->min_time_step_before_shutting_problematic_wells_ ==
179}
180
181template<class TypeTag>
182void
185{
186 registerEclTimeSteppingParameters<Scalar>();
188}
189
198template<class TypeTag>
199template <class Solver>
202step(const SimulatorTimer& simulator_timer,
203 Solver& solver,
204 const bool is_event,
205 const TuningUpdateCallback& tuning_updater)
206{
207 SubStepper<Solver> sub_stepper{
208 *this, simulator_timer, solver, is_event, tuning_updater,
209 };
210 return sub_stepper.run();
211}
212
213template<class TypeTag>
214template<class Serializer>
215void
217serializeOp(Serializer& serializer)
218{
219 serializer(this->time_step_control_type_);
220 switch (this->time_step_control_type_) {
222 allocAndSerialize<HardcodedTimeStepControl>(serializer);
223 break;
225 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
226 break;
228 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
229 break;
231 allocAndSerialize<PIDTimeStepControl>(serializer);
232 break;
234 allocAndSerialize<General3rdOrderController>(serializer);
235 break;
236 }
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_);
251}
252
253template<class TypeTag>
256report()
257{
258 return report_;
259}
260
261template<class TypeTag>
265{
266 return serializationTestObject_<HardcodedTimeStepControl>();
267}
268
269template<class TypeTag>
273{
274 return serializationTestObject_<PIDTimeStepControl>();
275}
276
277template<class TypeTag>
281{
282 return serializationTestObject_<PIDAndIterationCountTimeStepControl>();
283}
284
285template<class TypeTag>
289{
290 return serializationTestObject_<SimpleIterationCountTimeStepControl>();
291}
292
293template<class TypeTag>
297{
298 return serializationTestObject_<General3rdOrderController>();
299}
300
301
302template<class TypeTag>
303void
305setSuggestedNextStep(const double x)
306{
307 this->suggested_next_timestep_ = x;
308}
309
310template<class TypeTag>
311double
313suggestedNextStep() const
314{
315 return this->suggested_next_timestep_;
316}
317
318template<class TypeTag>
321timeStepControl() const
322{
323 return *this->time_step_control_;
324}
325
326
327template<class TypeTag>
328void
330updateNEXTSTEP(double max_next_tstep)
331{
332 // \Note Only update next suggested step if TSINIT was explicitly
333 // set in TUNING or NEXTSTEP is active.
334 if (max_next_tstep > 0) {
335 this->suggested_next_timestep_ = max_next_tstep;
336 }
337}
338
339template<class TypeTag>
340void
342updateTUNING(double max_next_tstep, const Tuning& tuning)
343{
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;
350}
351
352/*********************************************
353 * Private methods of AdaptiveTimeStepping
354 * ******************************************/
355
356template<class TypeTag>
357template<class T, class Serializer>
358void
360allocAndSerialize(Serializer& serializer)
361{
362 if (!serializer.isSerializing()) {
363 this->time_step_control_ = std::make_unique<T>();
364 }
365 serializer(*static_cast<T*>(this->time_step_control_.get()));
366}
367
368template<class TypeTag>
369template<class T>
370bool
371AdaptiveTimeStepping<TypeTag>::
372castAndComp(const AdaptiveTimeStepping<TypeTag>& Rhs) const
373{
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());
376 return *lhs == *rhs;
377}
378
379template<class TypeTag>
380void
381AdaptiveTimeStepping<TypeTag>::
382maybeModifySuggestedTimeStepAtBeginningOfReportStep_(const double original_time_step,
383 bool is_event)
384{
385 // init last time step as a fraction of the given time step
386 if (this->suggested_next_timestep_ < 0) {
387 this->suggested_next_timestep_ = this->restart_factor_ * original_time_step;
388 }
389
390 if (this->full_timestep_initially_) {
391 this->suggested_next_timestep_ = original_time_step;
392 }
393
394 // use seperate time step after event
395 if (is_event && this->timestep_after_event_ > 0) {
396 this->suggested_next_timestep_ = this->timestep_after_event_;
397 }
398}
399
400template<class TypeTag>
401template<class Controller>
402AdaptiveTimeStepping<TypeTag>
403AdaptiveTimeStepping<TypeTag>::
404serializationTestObject_()
405{
406 AdaptiveTimeStepping<TypeTag> result;
407
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());
424
425 return result;
426}
427
428/*********************************************
429 * Protected methods of AdaptiveTimeStepping
430 * ******************************************/
431
432template<class TypeTag>
434init_(const UnitSystem& unitSystem)
435{
436 std::tie(time_step_control_type_,
437 time_step_control_,
438 use_newton_iteration_) = detail::createController(unitSystem);
439 // make sure growth factor is something reasonable
440 if (this->growth_factor_ < 1.0) {
441 OPM_THROW(std::runtime_error,
442 "Growth factor cannot be less than 1.");
443 }
444}
445
446
447
448/************************************************
449 * Private class SubStepper public methods
450 ************************************************/
451
452template<class TypeTag>
453template<class Solver>
455SubStepper(AdaptiveTimeStepping<TypeTag>& adaptive_time_stepping,
456 const SimulatorTimer& simulator_timer,
457 Solver& solver,
458 const bool is_event,
459 const TuningUpdateCallback& tuning_updater)
460 : adaptive_time_stepping_{adaptive_time_stepping}
461 , simulator_timer_{simulator_timer}
462 , solver_{solver}
463 , is_event_{is_event}
464 , tuning_updater_{tuning_updater}
465{
466}
467
468template<class TypeTag>
469template<class Solver>
470AdaptiveTimeStepping<TypeTag>&
471AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
472getAdaptiveTimerStepper()
473{
474 return adaptive_time_stepping_;
475}
476
477template<class TypeTag>
478template<class Solver>
479SimulatorReport
480AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
481run()
482{
483#ifdef RESERVOIR_COUPLING_ENABLED
484 if (isReservoirCouplingSlave_() && reservoirCouplingSlave_().activated()) {
485 return runStepReservoirCouplingSlave_();
486 }
487 else if (isReservoirCouplingMaster_() && reservoirCouplingMaster_().activated()) {
488 return runStepReservoirCouplingMaster_();
489 }
490 else {
491 return runStepOriginal_();
492 }
493#else
494 return runStepOriginal_();
495#endif
496}
497
498/************************************************
499 * Private class SubStepper private methods
500 ************************************************/
501
502#ifdef RESERVOIR_COUPLING_ENABLED
503template<class TypeTag>
504template<class Solver>
505bool
506AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
507isReservoirCouplingMaster_() const
508{
509 return this->solver_.model().simulator().reservoirCouplingMaster() != nullptr;
510}
511
512template<class TypeTag>
513template<class Solver>
514bool
515AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
516isReservoirCouplingSlave_() const
517{
518 return this->solver_.model().simulator().reservoirCouplingSlave() != nullptr;
519}
520#endif
521
522template<class TypeTag>
523template<class Solver>
524void
525AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
526maybeModifySuggestedTimeStepAtBeginningOfReportStep_(const double original_time_step)
527{
528 this->adaptive_time_stepping_.maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
529 original_time_step, this->is_event_
530 );
531}
532
533// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
534// It has to be called for each substep since TUNING might have been changed for next sub step due
535// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
536template<class TypeTag>
537template<class Solver>
538bool
539AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
540maybeUpdateTuning_(double elapsed, double dt, int sub_step_number) const
541{
542 return this->tuning_updater_(elapsed, dt, sub_step_number);
543}
544
545template<class TypeTag>
546template<class Solver>
547double
548AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
549maxTimeStep_() const
550{
551 return this->adaptive_time_stepping_.max_time_step_;
552}
553
554template <class TypeTag>
555template <class Solver>
556SimulatorReport
557AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
558runStepOriginal_()
559{
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);
565
566 AdaptiveSimulatorTimer substep_timer{
567 this->simulator_timer_.startDateTime(),
568 original_time_step,
569 elapsed,
570 suggestedNextTimestep_(),
571 report_step,
572 maxTimeStep_()
573 };
574 SubStepIteration<Solver> substepIteration{*this, substep_timer, original_time_step, /*final_step=*/true};
575 return substepIteration.run();
576}
577
578#ifdef RESERVOIR_COUPLING_ENABLED
579template <class TypeTag>
580template <class Solver>
581ReservoirCouplingMaster<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
582AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
583reservoirCouplingMaster_()
584{
585 return *(this->solver_.model().simulator().reservoirCouplingMaster());
586}
587#endif
588
589#ifdef RESERVOIR_COUPLING_ENABLED
590template <class TypeTag>
591template <class Solver>
592ReservoirCouplingSlave<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
593AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
594reservoirCouplingSlave_()
595{
596 return *(this->solver_.model().simulator().reservoirCouplingSlave());
597}
598#endif
599
600#ifdef RESERVOIR_COUPLING_ENABLED
601// Description of the reservoir coupling master and slave substep loop
602// -------------------------------------------------------------------
603// The master and slave processes attempts to reach the end of the report step using a series of substeps
604// (also called timesteps). Each substep have an upper limit that is roughly determined by a combination
605// of the keywords TUNING (through the TSINIT, TSMAXZ values), NEXSTEP, WCYCLE, and the start of the
606// next report step (the last substep needs to coincide with this time). Note that NEXTSTEP can be
607// updated from an ACTIONX keyword. Although, this comment focuses on the maximum substep limit,
608// note that there is also a lower limit on the substep length. And the substep sizes will be adjusted
609// automatically (or retried) based on the convergence behavior of the solver and other criteria.
610//
611// When using reservoir coupling, the upper limit on each substep is further limited to: a) not overshoot
612// next report date of a slave reservoir, and b) to keep the flow rate of the slave groups within
613// certain limits. To determine this potential further limitation on the substep, the following procedure
614// is used at the beginning of each master substep:
615// - First, the slaves sends the master the date of their next report step
616// - The master receives the date of the next report step from the slaves
617// - If necessary, the master computes a modified substep length based on the received dates
618// TODO: explain how the substep is limited to keep the flow rate of the slave groups within certain
619// limits. Since this is not implemented yet, this part of the procedure is not explained here.
620//
621// Then, after the master has determined the substep length using the above procedure, it sends it
622// to the slaves. The slaves will now use the end of this substep as a fixed point (like a mini report
623// step), that they will try to reach using a single substep or multiple substeps. The end of the
624// substep received from the master will either conincide with the end of its current report step, or
625// it will end before (it cannot not end after since the master has already adjusted the step to not
626// exceed any report time in a slave). If the end of this substep is earlier in time than its next report
627// date, the slave will enter a loop and wait for the master to send it a new substep.
628// The loop is finished when the end of the received time step conincides with the end of its current
629// report step.
630
631template <class TypeTag>
632template <class Solver>
633SimulatorReport
634AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
635runStepReservoirCouplingMaster_()
636{
637 int iteration = 0;
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();
645 }
646 SimulatorReport report;
647 // The master needs to know which slaves have activated before it can start the substep loop
648 reservoirCouplingMaster_().maybeReceiveActivationHandshakeFromSlaves(current_time);
649 while (true) {
650 reservoirCouplingMaster_().receiveNextReportDateFromSlaves();
651 bool start_of_report_step = (iteration == 0);
652 if (start_of_report_step) {
653 reservoirCouplingMaster_().initStartOfReportStep(report_step_idx);
654 }
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, /*substep=*/0);
660 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
661 }
662 AdaptiveSimulatorTimer substep_timer{
663 this->simulator_timer_.startDateTime(),
664 /*stepLength=*/current_step_length,
665 /*elapsedTime=*/current_time,
666 /*timeStepEstimate=*/suggestedNextTimestep_(),
667 /*reportStep=*/this->simulator_timer_.reportStepNum(),
668 maxTimeStep_()
669 };
670 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
671 current_time + current_step_length, step_end_time
672 );
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;
677 if (final_step) {
678 break;
679 }
680 iteration++;
681 }
682 return report;
683}
684#endif
685
686#ifdef RESERVOIR_COUPLING_ENABLED
687template <class TypeTag>
688template <class Solver>
689SimulatorReport
690AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
691runStepReservoirCouplingSlave_()
692{
693 int iteration = 0;
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();
701 }
702 while (true) {
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, /*substep=*/0);
708 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
709 }
710 AdaptiveSimulatorTimer substep_timer{
711 this->simulator_timer_.startDateTime(),
712 /*step_length=*/timestep,
713 /*elapsed_time=*/current_time,
714 /*time_step_estimate=*/suggestedNextTimestep_(),
715 this->simulator_timer_.reportStepNum(),
716 maxTimeStep_()
717 };
718 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
719 current_time + timestep, step_end_time
720 );
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;
725 if (final_step) {
726 break;
727 }
728 iteration++;
729 }
730 return report;
731}
732
733#endif
734
735template <class TypeTag>
736template <class Solver>
737double
738AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
739suggestedNextTimestep_() const
740{
741 return this->adaptive_time_stepping_.suggestedNextStep();
742}
743
744
745
746/************************************************
747 * Private class SubStepIteration public methods
748 ************************************************/
749
750template<class TypeTag>
751template<class Solver>
752AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
753SubStepIteration(SubStepper<Solver>& substepper,
754 AdaptiveSimulatorTimer& substep_timer,
755 const double original_time_step,
756 bool final_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()}
762{
763}
764
765template <class TypeTag>
766template <class Solver>
767SimulatorReport
768AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
769run()
770{
771 auto& simulator = solver_().model().simulator();
772 auto& problem = simulator.problem();
773 // counter for solver restarts
774 int restarts = 0;
775 SimulatorReport report;
776
777 // sub step time loop
778 while (!this->substep_timer_.done()) {
779 // if we just chopped the timestep due to convergence i.e. restarts>0
780 // we dont want to update the next timestep based on Tuning
781 if (restarts == 0) {
782 maybeUpdateTuningAndTimeStep_();
783 }
784 const double dt = this->substep_timer_.currentStepLength();
785 if (timeStepVerbose_()) {
786 detail::logTimer(this->substep_timer_);
787 }
788
789 const auto substep_report = runSubStep_();
790
791 //Pass substep to eclwriter for summary output
792 problem.setSubStepReport(substep_report);
793 auto& full_report = adaptive_time_stepping_.report();
794 full_report += substep_report;
795 problem.setSimulationReport(full_report);
796
797 report += substep_report;
798
799 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
800 ++this->substep_timer_; // advance by current dt
801
802 const int iterations = getNumIterations_(substep_report);
803 auto dt_estimate = timeStepControlComputeEstimate_(
804 dt, iterations, this->substep_timer_);
805
806 assert(dt_estimate > 0);
807 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
808 restarts = 0; // solver converged, reset restarts counter
809
810 maybeReportSubStep_(substep_report);
811 if (this->final_step_ && this->substep_timer_.done()) {
812 // if the time step is done we do not need to write it as this will be done
813 // by the simulator anyway.
814 }
815 else {
816 report.success.output_write_time += writeOutput_();
817 }
818
819 // set new time step length
820 setTimeStep_(dt_estimate);
821
822 report.success.converged = this->substep_timer_.done();
823 this->substep_timer_.setLastStepFailed(false);
824 }
825 else { // in case of no convergence or time step tolerance test failure
826 this->substep_timer_.setLastStepFailed(true);
827 checkTimeStepMaxRestartLimit_(restarts);
828
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) { // added in case suggested time step is not a reduction
835 new_time_step = temp_time_step;
836 }
837 }
838 checkTimeStepMinLimit_(new_time_step);
839 bool wells_shut = false;
840 if (new_time_step > minTimeStepBeforeClosingWells_()) {
841 chopTimeStep_(new_time_step);
842 } else {
843 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
844 }
845 if (wells_shut) {
846 setTimeStep_(dt); // retry the old timestep
847 }
848 else {
849 restarts++; // only increase if no wells were shut
850 }
851 }
852 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
853 }
854 updateSuggestedNextStep_();
855 return report;
856}
857
858
859/************************************************
860 * Private class SubStepIteration private methods
861 ************************************************/
862
863
864template<class TypeTag>
865template<class Solver>
866bool
867AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
868checkContinueOnUnconvergedSolution_(double dt) const
869{
870 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
871 if (continue_on_uncoverged_solution && solverVerbose_()) {
872 // NOTE: This method is only called if the solver failed to converge.
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",
876 dt, minTimeStep_()
877 );
878 OpmLog::problem(msg);
879 }
880 return continue_on_uncoverged_solution;
881}
882
883template<class TypeTag>
884template<class Solver>
885void
886AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
887checkTimeStepMaxRestartLimit_(const int restarts) const
888{
889 // If we have restarted (i.e. cut the timestep) too
890 // many times, we have failed and throw an exception.
891 if (restarts >= solverRestartMax_()) {
892 const auto msg = fmt::format(
893 "Solver failed to converge after cutting timestep {} times.", restarts
894 );
895 if (solverVerbose_()) {
896 OpmLog::error(msg);
897 }
898 // Use throw directly to prevent file and line
899 throw TimeSteppingBreakdown{msg};
900 }
901}
902
903template<class TypeTag>
904template<class Solver>
905void
906AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
907checkTimeStepMinLimit_(const double new_time_step) const
908{
909 using Meas = UnitSystem::measure;
910 // If we have restarted (i.e. cut the timestep) too
911 // much, we have failed and throw an exception.
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();
916 msg += fmt::format(
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)
920 );
921 }
922 else {
923 msg += fmt::format(
924 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
925 minTimeStep_() / 86400.0
926 );
927 }
928 if (solverVerbose_()) {
929 OpmLog::error(msg);
930 }
931 // Use throw directly to prevent file and line
932 throw TimeSteppingBreakdown{msg};
933 }
934}
935
936template<class TypeTag>
937template<class Solver>
938void
939AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
940chopTimeStep_(const double new_time_step)
941{
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);
948 }
949}
950
951template<class TypeTag>
952template<class Solver>
953bool
954AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
955chopTimeStepOrCloseFailingWells_(const double new_time_step)
956{
957 bool wells_shut = false;
958 // We are below the threshold, and will check if there are any
959 // wells that fails repeatedly (that means that it fails in the last three steps)
960 // we should close rather than chopping again.
961 // If we already have chopped the timestep two times that is
962 // new_time_step < minTimeStepBeforeClosingWells_()*restartFactor_()*restartFactor_()
963 // We also shut wells that fails only on this step.
964 const bool requireRepeatedFailures =
965 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
966 const std::set<std::string> failing_wells =
967 detail::consistentlyFailingWells(solver_().model().stepReports(), requireRepeatedFailures);
968
969 if (failing_wells.empty()) {
970 // Found no wells to close, chop the timestep
971 chopTimeStep_(new_time_step);
972 } else {
973 // Close all consistently failing wells that are not under group control
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(),
979 /*dont_shut_grup_wells =*/ true);
980 if (was_shut) {
981 shut_wells.push_back(well);
982 }
983 }
984 // If no wells are closed we also try to shut wells under group control
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(),
990 /*dont_shut_grup_wells =*/ false);
991 if (was_shut) {
992 shut_wells.push_back(well);
993 }
994 }
995 }
996 // If still no wells are closed we must fall back to chopping again
997 if (shut_wells.empty()) {
998 chopTimeStep_(new_time_step);
999 } else {
1000 wells_shut = true;
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);
1007 }
1008 }
1009 }
1010 return wells_shut;
1011}
1012
1013template<class TypeTag>
1014template<class Solver>
1015boost::posix_time::ptime
1016AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1017currentDateTime_() const
1018{
1019 return simulatorTimer_().currentDateTime();
1020}
1021
1022template<class TypeTag>
1023template<class Solver>
1024int
1025AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1026getNumIterations_(const SimulatorReportSingle &substep_report) const
1027{
1028 if (useNewtonIteration_()) {
1029 return substep_report.total_newton_iterations;
1030 }
1031 else {
1032 return substep_report.total_linear_iterations;
1033 }
1034}
1035
1036template<class TypeTag>
1037template<class Solver>
1038double
1039AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1040growthFactor_() const
1041{
1042 return this->adaptive_time_stepping_.growth_factor_;
1043}
1044
1045template<class TypeTag>
1046template<class Solver>
1047bool
1048AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1049ignoreConvergenceFailure_() const
1050{
1051 return adaptive_time_stepping_.ignore_convergence_failure_;
1052}
1053
1054template<class TypeTag>
1055template<class Solver>
1056double
1057AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1058maxGrowth_() const
1059{
1060 return this->adaptive_time_stepping_.max_growth_;
1061}
1062
1063template<class TypeTag>
1064template<class Solver>
1065void
1066AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1067maybeReportSubStep_(SimulatorReportSingle substep_report) const
1068{
1069 if (timeStepVerbose_()) {
1070 std::ostringstream ss;
1071 substep_report.reportStep(ss);
1072 OpmLog::info(ss.str());
1073 }
1074}
1075
1076template<class TypeTag>
1077template<class Solver>
1078double
1079AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1080maybeRestrictTimeStepGrowth_(const double dt, double dt_estimate, const int restarts) const
1081{
1082 // limit the growth of the timestep size by the growth factor
1083 dt_estimate = std::min(dt_estimate, double(maxGrowth_() * dt));
1084 assert(dt_estimate > 0);
1085 // further restrict time step size growth after convergence problems
1086 if (restarts > 0) {
1087 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1088 }
1089
1090 return dt_estimate;
1091}
1092
1093// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
1094// It has to be called for each substep since TUNING might have been changed for next sub step due
1095// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
1096template<class TypeTag>
1097template<class Solver>
1098void
1099AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1100maybeUpdateTuningAndTimeStep_()
1101{
1102 // TODO: This function is currently only called if NEXTSTEP is activated from ACTIONX or
1103 // if the WCYCLE keyword needs to modify the current timestep. So this method should rather
1104 // be named maybeUpdateTimeStep_() or similar, since it should not update the tuning. However,
1105 // the current definition of the maybeUpdateTuning_() callback is actually calling
1106 // adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning) which is updating the tuning
1107 // see SimulatorFullyImplicitBlackoil::runStep() for more details.
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()))
1112 {
1113 // Either NEXTSTEP and WCYCLE wants to change the current time step, but they cannot
1114 // change the current time step directly. Instead, they change the suggested next time step
1115 // by calling updateNEXTSTEP() via the maybeUpdateTuning() callback. We now need to update
1116 // the current time step to the new suggested time step and reset the suggested time step
1117 // to the old value.
1118 setTimeStep_(suggestedNextTimestep_());
1119 setSuggestedNextStep_(old_value);
1120 }
1121}
1122
1123template<class TypeTag>
1124template<class Solver>
1125double
1126AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1127minTimeStepBeforeClosingWells_() const
1128{
1129 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1130}
1131
1132template<class TypeTag>
1133template<class Solver>
1134double
1135AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1136minTimeStep_() const
1137{
1138 return this->adaptive_time_stepping_.min_time_step_;
1139}
1140
1141template<class TypeTag>
1142template<class Solver>
1143double
1144AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1145restartFactor_() const
1146{
1147 return this->adaptive_time_stepping_.restart_factor_;
1148}
1149
1150template<class TypeTag>
1151template<class Solver>
1152SimulatorReportSingle
1153AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1154runSubStep_()
1155{
1156 SimulatorReportSingle substep_report;
1157
1158 auto handleFailure = [this, &substep_report]
1159 (const std::string& failure_reason, const std::exception& e, bool log_exception = true)
1160 {
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());
1165 }
1166 };
1167
1168 try {
1169 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1170 if (solverVerbose_()) {
1171 // report number of linear iterations
1172 OpmLog::debug("Overall linear iterations used: "
1173 + std::to_string(substep_report.total_linear_iterations));
1174 }
1175 }
1176 catch (const TooManyIterations& e) {
1177 handleFailure("Solver convergence failure - Iteration limit reached", e);
1178 }
1179 catch (const TimeSteppingBreakdown& e) {
1180 handleFailure(e.what(), e);
1181 }
1182 catch (const ConvergenceMonitorFailure& e) {
1183 handleFailure("Convergence monitor failure", e, /*log_exception=*/false);
1184 }
1185 catch (const LinearSolverProblem& e) {
1186 handleFailure("Linear solver convergence failure", e);
1187 }
1188 catch (const NumericalProblem& e) {
1189 handleFailure("Solver convergence failure - Numerical problem encountered", e);
1190 }
1191 catch (const std::runtime_error& e) {
1192 handleFailure("Runtime error encountered", e);
1193 }
1194 catch (const Dune::ISTLError& e) {
1195 handleFailure("ISTL error - Time step too large", e);
1196 }
1197 catch (const Dune::MatrixBlockError& e) {
1198 handleFailure("Matrix block error", e);
1199 }
1200
1201 return substep_report;
1202}
1203
1204template<class TypeTag>
1205template<class Solver>
1206void
1207AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1208setTimeStep_(double dt_estimate)
1209{
1210 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1211}
1212
1213template<class TypeTag>
1214template<class Solver>
1215Solver&
1216AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1217solver_() const
1218{
1219 return this->substepper_.solver_;
1220}
1221
1222
1223template<class TypeTag>
1224template<class Solver>
1225int
1226AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1227solverRestartMax_() const
1228{
1229 return this->adaptive_time_stepping_.solver_restart_max_;
1230}
1231
1232template<class TypeTag>
1233template<class Solver>
1234void
1235AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1236setSuggestedNextStep_(double step)
1237{
1238 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1239}
1240
1241template <class TypeTag>
1242template <class Solver>
1243const SimulatorTimer&
1244AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1245simulatorTimer_() const
1246{
1247 return this->substepper_.simulator_timer_;
1248}
1249
1250template <class TypeTag>
1251template <class Solver>
1252bool
1253AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1254solverVerbose_() const
1255{
1256 return this->adaptive_time_stepping_.solver_verbose_;
1257}
1258
1259template<class TypeTag>
1260template<class Solver>
1261boost::posix_time::ptime
1262AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1263startDateTime_() const
1264{
1265 return simulatorTimer_().startDateTime();
1266}
1267
1268template <class TypeTag>
1269template <class Solver>
1270double
1271AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1272suggestedNextTimestep_() const
1273{
1274 return this->adaptive_time_stepping_.suggestedNextStep();
1275}
1276
1277template <class TypeTag>
1278template <class Solver>
1279double
1280AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1281timeStepControlComputeEstimate_(const double dt, const int iterations,
1282 const AdaptiveSimulatorTimer& substepTimer) const
1283{
1284 // create object to compute the time error, simply forwards the call to the model
1285 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1286 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1287 dt, iterations, relative_change, substepTimer);
1288}
1289
1290template <class TypeTag>
1291template <class Solver>
1292bool
1293AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1294timeStepVerbose_() const
1295{
1296 return this->adaptive_time_stepping_.timestep_verbose_;
1297}
1298
1299// The suggested time step is the stepsize that will be used as a first try for
1300// the next sub step. It is updated at the end of each substep. It can also be
1301// updated by the TUNING or NEXTSTEP keywords at the beginning of each report step or
1302// at the beginning of each substep by the ACTIONX keyword (via NEXTSTEP), this is
1303// done by the maybeUpdateTuning_() method which is called at the beginning of each substep
1304// (and the begginning of each report step). Note that the WCYCLE keyword can also update the
1305// suggested time step via the maybeUpdateTuning_() method.
1306template <class TypeTag>
1307template <class Solver>
1308void
1309AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1310updateSuggestedNextStep_()
1311{
1312 auto suggested_next_step = this->substep_timer_.currentStepLength();
1313 if (! std::isfinite(suggested_next_step)) { // check for NaN
1314 suggested_next_step = this->original_time_step_;
1315 }
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());
1322 }
1323 setSuggestedNextStep_(suggested_next_step);
1324}
1325
1326template <class TypeTag>
1327template <class Solver>
1328bool
1329AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1330useNewtonIteration_() const
1331{
1332 return this->adaptive_time_stepping_.use_newton_iteration_;
1333}
1334
1335template <class TypeTag>
1336template <class Solver>
1337double
1338AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1339writeOutput_() const
1340{
1341 time::StopWatch perf_timer;
1342 perf_timer.start();
1343 auto& problem = solver_().model().simulator().problem();
1344 problem.writeOutput(true);
1345 return perf_timer.secsSinceStart();
1346}
1347
1348/************************************************
1349 * Private class SolutionTimeErrorSolverWrapper
1350 * **********************************************/
1351
1352template<class TypeTag>
1353template<class Solver>
1354AdaptiveTimeStepping<TypeTag>::
1355SolutionTimeErrorSolverWrapper<Solver>::
1356SolutionTimeErrorSolverWrapper(const Solver& solver)
1357 : solver_{solver}
1358{}
1359
1360template<class TypeTag>
1361template<class Solver>
1362double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange() const
1363{
1364 // returns: || u^n+1 - u^n || / || u^n+1 ||
1365 return solver_.model().relativeChange();
1366}
1367
1368} // namespace Opm
1369
1370#endif // OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
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
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