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&
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&
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 SimulatorReport report;
643 // The master needs to know which slaves have activated before it can start the substep loop
644 reservoirCouplingMaster_().maybeReceiveActivationHandshakeFromSlaves(current_time);
645 while (true) {
646 reservoirCouplingMaster_().receiveNextReportDateFromSlaves();
647 if (iteration == 0) {
648 maybeUpdateTuning_(current_time, current_step_length, /*substep=*/0);
649 }
650 current_step_length = reservoirCouplingMaster_().maybeChopSubStep(
651 current_step_length, current_time);
652 reservoirCouplingMaster_().sendNextTimeStepToSlaves(current_step_length);
653 if (iteration == 0) {
654 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
655 }
656 AdaptiveSimulatorTimer substep_timer{
657 this->simulator_timer_.startDateTime(),
658 /*stepLength=*/current_step_length,
659 /*elapsedTime=*/current_time,
660 /*timeStepEstimate=*/suggestedNextTimestep_(),
661 /*reportStep=*/this->simulator_timer_.reportStepNum(),
662 maxTimeStep_()
663 };
664 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
665 current_time + current_step_length, step_end_time
666 );
667 SubStepIteration<Solver> substepIteration{*this, substep_timer, current_step_length, final_step};
668 const auto sub_steps_report = substepIteration.run();
669 report += sub_steps_report;
670 current_time += current_step_length;
671 if (final_step) {
672 break;
673 }
674 iteration++;
675 }
676 return report;
677}
678#endif
679
680#ifdef RESERVOIR_COUPLING_ENABLED
681template <class TypeTag>
682template <class Solver>
683SimulatorReport
684AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
685runStepReservoirCouplingSlave_()
686{
687 int iteration = 0;
688 const double original_time_step = this->simulator_timer_.currentStepLength();
689 double current_time{this->simulator_timer_.simulationTimeElapsed()};
690 double step_end_time = current_time + original_time_step;
691 SimulatorReport report;
692 while (true) {
693 reservoirCouplingSlave_().sendNextReportDateToMasterProcess();
694 const auto timestep = reservoirCouplingSlave_().receiveNextTimeStepFromMaster();
695 if (iteration == 0) {
696 maybeUpdateTuning_(current_time, original_time_step, /*substep=*/0);
697 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
698 }
699 AdaptiveSimulatorTimer substep_timer{
700 this->simulator_timer_.startDateTime(),
701 /*step_length=*/timestep,
702 /*elapsed_time=*/current_time,
703 /*time_step_estimate=*/suggestedNextTimestep_(),
704 this->simulator_timer_.reportStepNum(),
705 maxTimeStep_()
706 };
707 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
708 current_time + timestep, step_end_time
709 );
710 SubStepIteration<Solver> substepIteration{*this, substep_timer, timestep, final_step};
711 const auto sub_steps_report = substepIteration.run();
712 report += sub_steps_report;
713 current_time += timestep;
714 if (final_step) {
715 break;
716 }
717 iteration++;
718 }
719 return report;
720}
721
722#endif
723
724template <class TypeTag>
725template <class Solver>
726double
727AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
728suggestedNextTimestep_() const
729{
730 return this->adaptive_time_stepping_.suggestedNextStep();
731}
732
733
734
735/************************************************
736 * Private class SubStepIteration public methods
737 ************************************************/
738
739template<class TypeTag>
740template<class Solver>
741AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
742SubStepIteration(SubStepper<Solver>& substepper,
743 AdaptiveSimulatorTimer& substep_timer,
744 const double original_time_step,
745 bool final_step)
746 : substepper_{substepper}
747 , substep_timer_{substep_timer}
748 , original_time_step_{original_time_step}
749 , final_step_{final_step}
750 , adaptive_time_stepping_{substepper.getAdaptiveTimerStepper()}
751{
752}
753
754template <class TypeTag>
755template <class Solver>
756SimulatorReport
757AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
758run()
759{
760 auto& simulator = solver_().model().simulator();
761 auto& problem = simulator.problem();
762 // counter for solver restarts
763 int restarts = 0;
764 SimulatorReport report;
765
766 // sub step time loop
767 while (!this->substep_timer_.done()) {
768 // if we just chopped the timestep due to convergence i.e. restarts>0
769 // we dont want to update the next timestep based on Tuning
770 if (restarts == 0) {
771 maybeUpdateTuningAndTimeStep_();
772 }
773 const double dt = this->substep_timer_.currentStepLength();
774 if (timeStepVerbose_()) {
775 detail::logTimer(this->substep_timer_);
776 }
777
778 const auto substep_report = runSubStep_();
779
780 //Pass substep to eclwriter for summary output
781 problem.setSubStepReport(substep_report);
782 auto& full_report = adaptive_time_stepping_.report();
783 full_report += substep_report;
784 problem.setSimulationReport(full_report);
785
786 report += substep_report;
787
788 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
789 ++this->substep_timer_; // advance by current dt
790
791 const int iterations = getNumIterations_(substep_report);
792 auto dt_estimate = timeStepControlComputeEstimate_(
793 dt, iterations, this->substep_timer_);
794
795 assert(dt_estimate > 0);
796 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
797 restarts = 0; // solver converged, reset restarts counter
798
799 maybeReportSubStep_(substep_report);
800 if (this->final_step_ && this->substep_timer_.done()) {
801 // if the time step is done we do not need to write it as this will be done
802 // by the simulator anyway.
803 }
804 else {
805 report.success.output_write_time += writeOutput_();
806 }
807
808 // set new time step length
809 setTimeStep_(dt_estimate);
810
811 report.success.converged = this->substep_timer_.done();
812 this->substep_timer_.setLastStepFailed(false);
813 }
814 else { // in case of no convergence or time step tolerance test failure
815 this->substep_timer_.setLastStepFailed(true);
816 checkTimeStepMaxRestartLimit_(restarts);
817
818 double new_time_step = restartFactor_() * dt;
819 if (substep_report.time_step_rejected) {
820 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>();
821 const double safetyFactor = Parameters::Get<Parameters::TimeStepControlSafetyFactor>();
822 const double temp_time_step = std::sqrt(safetyFactor * tol / solver_().model().relativeChange()) * dt;
823 if (temp_time_step < dt) { // added in case suggested time step is not a reduction
824 new_time_step = temp_time_step;
825 }
826 }
827 checkTimeStepMinLimit_(new_time_step);
828 bool wells_shut = false;
829 if (new_time_step > minTimeStepBeforeClosingWells_()) {
830 chopTimeStep_(new_time_step);
831 } else {
832 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
833 }
834 if (wells_shut) {
835 setTimeStep_(dt); // retry the old timestep
836 }
837 else {
838 restarts++; // only increase if no wells were shut
839 }
840 }
841 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
842 }
843 updateSuggestedNextStep_();
844 return report;
845}
846
847
848/************************************************
849 * Private class SubStepIteration private methods
850 ************************************************/
851
852
853template<class TypeTag>
854template<class Solver>
855bool
856AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
857checkContinueOnUnconvergedSolution_(double dt) const
858{
859 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
860 if (continue_on_uncoverged_solution && solverVerbose_()) {
861 // NOTE: This method is only called if the solver failed to converge.
862 const auto msg = fmt::format(
863 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
864 "which is the minimum threshold given by option --solver-min-time-step\n",
865 dt, minTimeStep_()
866 );
867 OpmLog::problem(msg);
868 }
869 return continue_on_uncoverged_solution;
870}
871
872template<class TypeTag>
873template<class Solver>
874void
875AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
876checkTimeStepMaxRestartLimit_(const int restarts) const
877{
878 // If we have restarted (i.e. cut the timestep) too
879 // many times, we have failed and throw an exception.
880 if (restarts >= solverRestartMax_()) {
881 const auto msg = fmt::format(
882 "Solver failed to converge after cutting timestep {} times.", restarts
883 );
884 if (solverVerbose_()) {
885 OpmLog::error(msg);
886 }
887 // Use throw directly to prevent file and line
888 throw TimeSteppingBreakdown{msg};
889 }
890}
891
892template<class TypeTag>
893template<class Solver>
894void
895AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
896checkTimeStepMinLimit_(const double new_time_step) const
897{
898 using Meas = UnitSystem::measure;
899 // If we have restarted (i.e. cut the timestep) too
900 // much, we have failed and throw an exception.
901 if (new_time_step < minTimeStep_()) {
902 auto msg = fmt::format("Solver failed to converge after cutting timestep to ");
903 if (Parameters::Get<Parameters::EnableTuning>()) {
904 const UnitSystem& unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
905 msg += fmt::format(
906 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
907 unit_system.from_si(Meas::time, minTimeStep_()),
908 unit_system.name(Meas::time)
909 );
910 }
911 else {
912 msg += fmt::format(
913 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
914 minTimeStep_() / 86400.0
915 );
916 }
917 if (solverVerbose_()) {
918 OpmLog::error(msg);
919 }
920 // Use throw directly to prevent file and line
921 throw TimeSteppingBreakdown{msg};
922 }
923}
924
925template<class TypeTag>
926template<class Solver>
927void
928AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
929chopTimeStep_(const double new_time_step)
930{
931 setTimeStep_(new_time_step);
932 if (solverVerbose_()) {
933 const auto msg = fmt::format("{}\nTimestep chopped to {} days\n",
934 this->cause_of_failure_,
935 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
936 OpmLog::problem(msg);
937 }
938}
939
940template<class TypeTag>
941template<class Solver>
942bool
943AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
944chopTimeStepOrCloseFailingWells_(const double new_time_step)
945{
946 bool wells_shut = false;
947 // We are below the threshold, and will check if there are any
948 // wells that fails repeatedly (that means that it fails in the last three steps)
949 // we should close rather than chopping again.
950 // If we already have chopped the timestep two times that is
951 // new_time_step < minTimeStepBeforeClosingWells_()*restartFactor_()*restartFactor_()
952 // We also shut wells that fails only on this step.
953 const bool requireRepeatedFailures =
954 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
955 const std::set<std::string> failing_wells =
956 detail::consistentlyFailingWells(solver_().model().stepReports(), requireRepeatedFailures);
957
958 if (failing_wells.empty()) {
959 // Found no wells to close, chop the timestep
960 chopTimeStep_(new_time_step);
961 } else {
962 // Close all consistently failing wells that are not under group control
963 std::vector<std::string> shut_wells;
964 for (const auto& well : failing_wells) {
965 const bool was_shut =
966 solver_().model().wellModel().forceShutWellByName(well,
967 this->substep_timer_.simulationTimeElapsed(),
968 /*dont_shut_grup_wells =*/ true);
969 if (was_shut) {
970 shut_wells.push_back(well);
971 }
972 }
973 // If no wells are closed we also try to shut wells under group control
974 if (shut_wells.empty()) {
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 =*/ false);
980 if (was_shut) {
981 shut_wells.push_back(well);
982 }
983 }
984 }
985 // If still no wells are closed we must fall back to chopping again
986 if (shut_wells.empty()) {
987 chopTimeStep_(new_time_step);
988 } else {
989 wells_shut = true;
990 if (solverVerbose_()) {
991 const std::string msg =
992 fmt::format("\nProblematic well(s) were shut: {}"
993 "(retrying timestep)\n",
994 fmt::join(shut_wells, " "));
995 OpmLog::problem(msg);
996 }
997 }
998 }
999 return wells_shut;
1000}
1001
1002template<class TypeTag>
1003template<class Solver>
1004boost::posix_time::ptime
1005AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1006currentDateTime_() const
1007{
1008 return simulatorTimer_().currentDateTime();
1009}
1010
1011template<class TypeTag>
1012template<class Solver>
1013int
1014AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1015getNumIterations_(const SimulatorReportSingle &substep_report) const
1016{
1017 if (useNewtonIteration_()) {
1018 return substep_report.total_newton_iterations;
1019 }
1020 else {
1021 return substep_report.total_linear_iterations;
1022 }
1023}
1024
1025template<class TypeTag>
1026template<class Solver>
1027double
1028AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1029growthFactor_() const
1030{
1031 return this->adaptive_time_stepping_.growth_factor_;
1032}
1033
1034template<class TypeTag>
1035template<class Solver>
1036bool
1037AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1038ignoreConvergenceFailure_() const
1039{
1040 return adaptive_time_stepping_.ignore_convergence_failure_;
1041}
1042
1043template<class TypeTag>
1044template<class Solver>
1045double
1046AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1047maxGrowth_() const
1048{
1049 return this->adaptive_time_stepping_.max_growth_;
1050}
1051
1052template<class TypeTag>
1053template<class Solver>
1054void
1055AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1056maybeReportSubStep_(SimulatorReportSingle substep_report) const
1057{
1058 if (timeStepVerbose_()) {
1059 std::ostringstream ss;
1060 substep_report.reportStep(ss);
1061 OpmLog::info(ss.str());
1062 }
1063}
1064
1065template<class TypeTag>
1066template<class Solver>
1067double
1068AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1069maybeRestrictTimeStepGrowth_(const double dt, double dt_estimate, const int restarts) const
1070{
1071 // limit the growth of the timestep size by the growth factor
1072 dt_estimate = std::min(dt_estimate, double(maxGrowth_() * dt));
1073 assert(dt_estimate > 0);
1074 // further restrict time step size growth after convergence problems
1075 if (restarts > 0) {
1076 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1077 }
1078
1079 return dt_estimate;
1080}
1081
1082// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
1083// It has to be called for each substep since TUNING might have been changed for next sub step due
1084// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
1085template<class TypeTag>
1086template<class Solver>
1087void
1088AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1089maybeUpdateTuningAndTimeStep_()
1090{
1091 // TODO: This function is currently only called if NEXTSTEP is activated from ACTIONX or
1092 // if the WCYCLE keyword needs to modify the current timestep. So this method should rather
1093 // be named maybeUpdateTimeStep_() or similar, since it should not update the tuning. However,
1094 // the current definition of the maybeUpdateTuning_() callback is actually calling
1095 // adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning) which is updating the tuning
1096 // see SimulatorFullyImplicitBlackoil::runStep() for more details.
1097 const auto old_value = suggestedNextTimestep_();
1098 if (this->substepper_.maybeUpdateTuning_(this->substep_timer_.simulationTimeElapsed(),
1099 this->substep_timer_.currentStepLength(),
1100 this->substep_timer_.currentStepNum()))
1101 {
1102 // Either NEXTSTEP and WCYCLE wants to change the current time step, but they cannot
1103 // change the current time step directly. Instead, they change the suggested next time step
1104 // by calling updateNEXTSTEP() via the maybeUpdateTuning() callback. We now need to update
1105 // the current time step to the new suggested time step and reset the suggested time step
1106 // to the old value.
1107 setTimeStep_(suggestedNextTimestep_());
1108 setSuggestedNextStep_(old_value);
1109 }
1110}
1111
1112template<class TypeTag>
1113template<class Solver>
1114double
1115AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1116minTimeStepBeforeClosingWells_() const
1117{
1118 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1119}
1120
1121template<class TypeTag>
1122template<class Solver>
1123double
1124AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1125minTimeStep_() const
1126{
1127 return this->adaptive_time_stepping_.min_time_step_;
1128}
1129
1130template<class TypeTag>
1131template<class Solver>
1132double
1133AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1134restartFactor_() const
1135{
1136 return this->adaptive_time_stepping_.restart_factor_;
1137}
1138
1139template<class TypeTag>
1140template<class Solver>
1141SimulatorReportSingle
1142AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1143runSubStep_()
1144{
1145 SimulatorReportSingle substep_report;
1146
1147 auto handleFailure = [this, &substep_report]
1148 (const std::string& failure_reason, const std::exception& e, bool log_exception = true)
1149 {
1150 substep_report = solver_().failureReport();
1151 this->cause_of_failure_ = failure_reason;
1152 if (log_exception && solverVerbose_()) {
1153 OpmLog::debug(std::string("Caught Exception: ") + e.what());
1154 }
1155 };
1156
1157 try {
1158 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1159 if (solverVerbose_()) {
1160 // report number of linear iterations
1161 OpmLog::debug("Overall linear iterations used: "
1162 + std::to_string(substep_report.total_linear_iterations));
1163 }
1164 }
1165 catch (const TooManyIterations& e) {
1166 handleFailure("Solver convergence failure - Iteration limit reached", e);
1167 }
1168 catch (const TimeSteppingBreakdown& e) {
1169 handleFailure(e.what(), e);
1170 }
1171 catch (const ConvergenceMonitorFailure& e) {
1172 handleFailure("Convergence monitor failure", e, /*log_exception=*/false);
1173 }
1174 catch (const LinearSolverProblem& e) {
1175 handleFailure("Linear solver convergence failure", e);
1176 }
1177 catch (const NumericalProblem& e) {
1178 handleFailure("Solver convergence failure - Numerical problem encountered", e);
1179 }
1180 catch (const std::runtime_error& e) {
1181 handleFailure("Runtime error encountered", e);
1182 }
1183 catch (const Dune::ISTLError& e) {
1184 handleFailure("ISTL error - Time step too large", e);
1185 }
1186 catch (const Dune::MatrixBlockError& e) {
1187 handleFailure("Matrix block error", e);
1188 }
1189
1190 return substep_report;
1191}
1192
1193template<class TypeTag>
1194template<class Solver>
1195void
1196AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1197setTimeStep_(double dt_estimate)
1198{
1199 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1200}
1201
1202template<class TypeTag>
1203template<class Solver>
1204Solver&
1205AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1206solver_() const
1207{
1208 return this->substepper_.solver_;
1209}
1210
1211
1212template<class TypeTag>
1213template<class Solver>
1214int
1215AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1216solverRestartMax_() const
1217{
1218 return this->adaptive_time_stepping_.solver_restart_max_;
1219}
1220
1221template<class TypeTag>
1222template<class Solver>
1223void
1224AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1225setSuggestedNextStep_(double step)
1226{
1227 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1228}
1229
1230template <class TypeTag>
1231template <class Solver>
1232const SimulatorTimer&
1233AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1234simulatorTimer_() const
1235{
1236 return this->substepper_.simulator_timer_;
1237}
1238
1239template <class TypeTag>
1240template <class Solver>
1241bool
1242AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1243solverVerbose_() const
1244{
1245 return this->adaptive_time_stepping_.solver_verbose_;
1246}
1247
1248template<class TypeTag>
1249template<class Solver>
1250boost::posix_time::ptime
1251AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1252startDateTime_() const
1253{
1254 return simulatorTimer_().startDateTime();
1255}
1256
1257template <class TypeTag>
1258template <class Solver>
1259double
1260AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1261suggestedNextTimestep_() const
1262{
1263 return this->adaptive_time_stepping_.suggestedNextStep();
1264}
1265
1266template <class TypeTag>
1267template <class Solver>
1268double
1269AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1270timeStepControlComputeEstimate_(const double dt, const int iterations,
1271 const AdaptiveSimulatorTimer& substepTimer) const
1272{
1273 // create object to compute the time error, simply forwards the call to the model
1274 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1275 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1276 dt, iterations, relative_change, substepTimer);
1277}
1278
1279template <class TypeTag>
1280template <class Solver>
1281bool
1282AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1283timeStepVerbose_() const
1284{
1285 return this->adaptive_time_stepping_.timestep_verbose_;
1286}
1287
1288// The suggested time step is the stepsize that will be used as a first try for
1289// the next sub step. It is updated at the end of each substep. It can also be
1290// updated by the TUNING or NEXTSTEP keywords at the beginning of each report step or
1291// at the beginning of each substep by the ACTIONX keyword (via NEXTSTEP), this is
1292// done by the maybeUpdateTuning_() method which is called at the beginning of each substep
1293// (and the begginning of each report step). Note that the WCYCLE keyword can also update the
1294// suggested time step via the maybeUpdateTuning_() method.
1295template <class TypeTag>
1296template <class Solver>
1297void
1298AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1299updateSuggestedNextStep_()
1300{
1301 auto suggested_next_step = this->substep_timer_.currentStepLength();
1302 if (! std::isfinite(suggested_next_step)) { // check for NaN
1303 suggested_next_step = this->original_time_step_;
1304 }
1305 if (timeStepVerbose_()) {
1306 std::ostringstream ss;
1307 this->substep_timer_.report(ss);
1308 ss << "Suggested next step size = "
1309 << unit::convert::to(suggested_next_step, unit::day) << " (days)" << std::endl;
1310 OpmLog::debug(ss.str());
1311 }
1312 setSuggestedNextStep_(suggested_next_step);
1313}
1314
1315template <class TypeTag>
1316template <class Solver>
1317bool
1318AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1319useNewtonIteration_() const
1320{
1321 return this->adaptive_time_stepping_.use_newton_iteration_;
1322}
1323
1324template <class TypeTag>
1325template <class Solver>
1326double
1327AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1328writeOutput_() const
1329{
1330 time::StopWatch perf_timer;
1331 perf_timer.start();
1332 auto& problem = solver_().model().simulator().problem();
1333 problem.writeOutput(true);
1334 return perf_timer.secsSinceStart();
1335}
1336
1337/************************************************
1338 * Private class SolutionTimeErrorSolverWrapper
1339 * **********************************************/
1340
1341template<class TypeTag>
1342template<class Solver>
1343AdaptiveTimeStepping<TypeTag>::
1344SolutionTimeErrorSolverWrapper<Solver>::
1345SolutionTimeErrorSolverWrapper(const Solver& solver)
1346 : solver_{solver}
1347{}
1348
1349template<class TypeTag>
1350template<class Solver>
1351double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange() const
1352{
1353 // returns: || u^n+1 - u^n || / || u^n+1 ||
1354 return solver_.model().relativeChange();
1355}
1356
1357} // namespace Opm
1358
1359#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
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
bool operator==(const AdaptiveTimeStepping< TypeTag > &rhs)
Definition: AdaptiveTimeStepping_impl.hpp:138
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: blackoilboundaryratevector.hh:39
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
This file provides the infrastructure to retrieve run-time parameters.
static bool compare_gt_or_eq(double a, double b)
Determines if a is greater than b within the specified tolerance.
Definition: SimulatorReport.hpp:122