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
100template<class TypeTag>
102AdaptiveTimeStepping(double max_next_tstep,
103 const Tuning& tuning,
104 const UnitSystem& unit_system,
105 const SimulatorReport& report,
106 const bool terminal_output
107)
108 : time_step_control_{}
109 , restart_factor_{tuning.TSFCNV}
110 , growth_factor_{tuning.TFDIFF}
111 , max_growth_{tuning.TSFMAX}
112 , max_time_step_{tuning.TSMAXZ} // 365.25
113 , min_time_step_{tuning.TSMINZ} // 0.1;
114 , ignore_convergence_failure_{true}
115 , solver_restart_max_{Parameters::Get<Parameters::SolverMaxRestarts>()} // 10
116 , solver_verbose_{Parameters::Get<Parameters::SolverVerbosity>() > 0 && terminal_output} // 2
117 , timestep_verbose_{Parameters::Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output} // 2
118 , suggested_next_timestep_{
119 max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>() * 24 * 60 * 60
120 : max_next_tstep} // 1.0
121 , full_timestep_initially_{Parameters::Get<Parameters::FullTimeStepInitially>()} // false
122 , timestep_after_event_{tuning.TMAXWC} // 1e30
123 , use_newton_iteration_{false}
124 , min_time_step_before_shutting_problematic_wells_{
125 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
126 , report_(report)
127{
128 init_(unit_system);
129}
130
131template<class TypeTag>
132bool
135{
136 if (this->time_step_control_type_ != rhs.time_step_control_type_ ||
137 (this->time_step_control_ && !rhs.time_step_control_) ||
138 (!this->time_step_control_ && rhs.time_step_control_)) {
139 return false;
140 }
141
142 bool result = false;
143 switch (this->time_step_control_type_) {
145 result = castAndComp<HardcodedTimeStepControl>(rhs);
146 break;
148 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
149 break;
151 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
152 break;
154 result = castAndComp<PIDTimeStepControl>(rhs);
155 break;
157 result = castAndComp<General3rdOrderController>(rhs);
158 break;
159 }
160
161 return result &&
162 this->restart_factor_ == rhs.restart_factor_ &&
163 this->growth_factor_ == rhs.growth_factor_ &&
164 this->max_growth_ == rhs.max_growth_ &&
165 this->max_time_step_ == rhs.max_time_step_ &&
166 this->min_time_step_ == rhs.min_time_step_ &&
167 this->ignore_convergence_failure_ == rhs.ignore_convergence_failure_ &&
168 this->solver_restart_max_== rhs.solver_restart_max_ &&
169 this->solver_verbose_ == rhs.solver_verbose_ &&
170 this->full_timestep_initially_ == rhs.full_timestep_initially_ &&
171 this->timestep_after_event_ == rhs.timestep_after_event_ &&
172 this->use_newton_iteration_ == rhs.use_newton_iteration_ &&
173 this->min_time_step_before_shutting_problematic_wells_ ==
175}
176
177template<class TypeTag>
178void
181{
182 registerEclTimeSteppingParameters<Scalar>();
184}
185
191template<class TypeTag>
192template <class Solver>
195step(const SimulatorTimer& simulator_timer,
196 Solver& solver,
197 const bool is_event,
198 const TuningUpdateCallback& tuning_updater)
199{
200 SubStepper<Solver> sub_stepper{
201 *this, simulator_timer, solver, is_event, tuning_updater,
202 };
203 return sub_stepper.run();
204}
205
206template<class TypeTag>
207template<class Serializer>
208void
210serializeOp(Serializer& serializer)
211{
212 serializer(this->time_step_control_type_);
213 switch (this->time_step_control_type_) {
215 allocAndSerialize<HardcodedTimeStepControl>(serializer);
216 break;
218 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
219 break;
221 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
222 break;
224 allocAndSerialize<PIDTimeStepControl>(serializer);
225 break;
227 allocAndSerialize<General3rdOrderController>(serializer);
228 break;
229 }
230 serializer(this->restart_factor_);
231 serializer(this->growth_factor_);
232 serializer(this->max_growth_);
233 serializer(this->max_time_step_);
234 serializer(this->min_time_step_);
235 serializer(this->ignore_convergence_failure_);
236 serializer(this->solver_restart_max_);
237 serializer(this->solver_verbose_);
238 serializer(this->timestep_verbose_);
239 serializer(this->suggested_next_timestep_);
240 serializer(this->full_timestep_initially_);
241 serializer(this->timestep_after_event_);
242 serializer(this->use_newton_iteration_);
243 serializer(this->min_time_step_before_shutting_problematic_wells_);
244}
245
246template<class TypeTag>
249report()
250{
251 return report_;
252}
253
254template<class TypeTag>
258{
259 return serializationTestObject_<HardcodedTimeStepControl>();
260}
261
262template<class TypeTag>
266{
267 return serializationTestObject_<PIDTimeStepControl>();
268}
269
270template<class TypeTag>
274{
275 return serializationTestObject_<PIDAndIterationCountTimeStepControl>();
276}
277
278template<class TypeTag>
282{
283 return serializationTestObject_<SimpleIterationCountTimeStepControl>();
284}
285
286template<class TypeTag>
290{
291 return serializationTestObject_<General3rdOrderController>();
292}
293
294
295template<class TypeTag>
296void
298setSuggestedNextStep(const double x)
299{
300 this->suggested_next_timestep_ = x;
301}
302
303template<class TypeTag>
304double
306suggestedNextStep() const
307{
308 return this->suggested_next_timestep_;
309}
310
311template<class TypeTag>
314timeStepControl() const
315{
316 return *this->time_step_control_;
317}
318
319
320template<class TypeTag>
321void
323updateNEXTSTEP(double max_next_tstep)
324{
325 // \Note Only update next suggested step if TSINIT was explicitly
326 // set in TUNING or NEXTSTEP is active.
327 if (max_next_tstep > 0) {
328 this->suggested_next_timestep_ = max_next_tstep;
329 }
330}
331
332template<class TypeTag>
333void
335updateTUNING(double max_next_tstep, const Tuning& tuning)
336{
337 this->restart_factor_ = tuning.TSFCNV;
338 this->growth_factor_ = tuning.TFDIFF;
339 this->max_growth_ = tuning.TSFMAX;
340 this->max_time_step_ = tuning.TSMAXZ;
341 updateNEXTSTEP(max_next_tstep);
342 this->timestep_after_event_ = tuning.TMAXWC;
343}
344
345/*********************************************
346 * Private methods of AdaptiveTimeStepping
347 * ******************************************/
348
349template<class TypeTag>
350template<class T, class Serializer>
351void
353allocAndSerialize(Serializer& serializer)
354{
355 if (!serializer.isSerializing()) {
356 this->time_step_control_ = std::make_unique<T>();
357 }
358 serializer(*static_cast<T*>(this->time_step_control_.get()));
359}
360
361template<class TypeTag>
362template<class T>
363bool
364AdaptiveTimeStepping<TypeTag>::
365castAndComp(const AdaptiveTimeStepping<TypeTag>& Rhs) const
366{
367 const T* lhs = static_cast<const T*>(this->time_step_control_.get());
368 const T* rhs = static_cast<const T*>(Rhs.time_step_control_.get());
369 return *lhs == *rhs;
370}
371
372template<class TypeTag>
373void
374AdaptiveTimeStepping<TypeTag>::
375maybeModifySuggestedTimeStepAtBeginningOfReportStep_(const double original_time_step,
376 bool is_event)
377{
378 // init last time step as a fraction of the given time step
379 if (this->suggested_next_timestep_ < 0) {
380 this->suggested_next_timestep_ = this->restart_factor_ * original_time_step;
381 }
382
383 if (this->full_timestep_initially_) {
384 this->suggested_next_timestep_ = original_time_step;
385 }
386
387 // use seperate time step after event
388 if (is_event && this->timestep_after_event_ > 0) {
389 this->suggested_next_timestep_ = this->timestep_after_event_;
390 }
391}
392
393template<class TypeTag>
394template<class Controller>
395AdaptiveTimeStepping<TypeTag>
396AdaptiveTimeStepping<TypeTag>::
397serializationTestObject_()
398{
399 AdaptiveTimeStepping<TypeTag> result;
400
401 result.restart_factor_ = 1.0;
402 result.growth_factor_ = 2.0;
403 result.max_growth_ = 3.0;
404 result.max_time_step_ = 4.0;
405 result.min_time_step_ = 5.0;
406 result.ignore_convergence_failure_ = true;
407 result.solver_restart_max_ = 6;
408 result.solver_verbose_ = true;
409 result.timestep_verbose_ = true;
410 result.suggested_next_timestep_ = 7.0;
411 result.full_timestep_initially_ = true;
412 result.use_newton_iteration_ = true;
413 result.min_time_step_before_shutting_problematic_wells_ = 9.0;
414 result.time_step_control_type_ = Controller::Type;
415 result.time_step_control_ =
416 std::make_unique<Controller>(Controller::serializationTestObject());
417
418 return result;
419}
420
421/*********************************************
422 * Protected methods of AdaptiveTimeStepping
423 * ******************************************/
424
425template<class TypeTag>
427init_(const UnitSystem& unitSystem)
428{
429 std::tie(time_step_control_type_,
430 time_step_control_,
431 use_newton_iteration_) = detail::createController(unitSystem);
432 // make sure growth factor is something reasonable
433 if (this->growth_factor_ < 1.0) {
434 OPM_THROW(std::runtime_error,
435 "Growth factor cannot be less than 1.");
436 }
437}
438
439
440
441/************************************************
442 * Private class SubStepper public methods
443 ************************************************/
444
445template<class TypeTag>
446template<class Solver>
448SubStepper(AdaptiveTimeStepping<TypeTag>& adaptive_time_stepping,
449 const SimulatorTimer& simulator_timer,
450 Solver& solver,
451 const bool is_event,
452 const TuningUpdateCallback& tuning_updater)
453 : adaptive_time_stepping_{adaptive_time_stepping}
454 , simulator_timer_{simulator_timer}
455 , solver_{solver}
456 , is_event_{is_event}
457 , tuning_updater_{tuning_updater}
458{
459}
460
461template<class TypeTag>
462template<class Solver>
463AdaptiveTimeStepping<TypeTag>&
464AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
465getAdaptiveTimerStepper()
466{
467 return adaptive_time_stepping_;
468}
469
470template<class TypeTag>
471template<class Solver>
472SimulatorReport
473AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
474run()
475{
476#ifdef RESERVOIR_COUPLING_ENABLED
477 if (isReservoirCouplingSlave_() && reservoirCouplingSlave_().activated()) {
478 return runStepReservoirCouplingSlave_();
479 }
480 else if (isReservoirCouplingMaster_() && reservoirCouplingMaster_().activated()) {
481 return runStepReservoirCouplingMaster_();
482 }
483 else {
484 return runStepOriginal_();
485 }
486#else
487 return runStepOriginal_();
488#endif
489}
490
491/************************************************
492 * Private class SubStepper private methods
493 ************************************************/
494
495#ifdef RESERVOIR_COUPLING_ENABLED
496template<class TypeTag>
497template<class Solver>
498bool
499AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
500isReservoirCouplingMaster_() const
501{
502 return this->solver_.model().simulator().reservoirCouplingMaster() != nullptr;
503}
504
505template<class TypeTag>
506template<class Solver>
507bool
508AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
509isReservoirCouplingSlave_() const
510{
511 return this->solver_.model().simulator().reservoirCouplingSlave() != nullptr;
512}
513#endif
514
515template<class TypeTag>
516template<class Solver>
517void
518AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
519maybeModifySuggestedTimeStepAtBeginningOfReportStep_(const double original_time_step)
520{
521 this->adaptive_time_stepping_.maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
522 original_time_step, this->is_event_
523 );
524}
525
526// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
527// It has to be called for each substep since TUNING might have been changed for next sub step due
528// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
529template<class TypeTag>
530template<class Solver>
531bool
532AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
533maybeUpdateTuning_(double elapsed, double dt, int sub_step_number) const
534{
535 return this->tuning_updater_(elapsed, dt, sub_step_number);
536}
537
538template<class TypeTag>
539template<class Solver>
540double
541AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
542maxTimeStep_() const
543{
544 return this->adaptive_time_stepping_.max_time_step_;
545}
546
547template <class TypeTag>
548template <class Solver>
549SimulatorReport
550AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
551runStepOriginal_()
552{
553 const auto elapsed = this->simulator_timer_.simulationTimeElapsed();
554 const auto original_time_step = this->simulator_timer_.currentStepLength();
555 const auto report_step = this->simulator_timer_.reportStepNum();
556 maybeUpdateTuning_(elapsed, original_time_step, report_step);
557 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(original_time_step);
558
559 AdaptiveSimulatorTimer substep_timer{
560 this->simulator_timer_.startDateTime(),
561 original_time_step,
562 elapsed,
563 suggestedNextTimestep_(),
564 report_step,
565 maxTimeStep_()
566 };
567 SubStepIteration<Solver> substepIteration{*this, substep_timer, original_time_step, /*final_step=*/true};
568 return substepIteration.run();
569}
570
571#ifdef RESERVOIR_COUPLING_ENABLED
572template <class TypeTag>
573template <class Solver>
574ReservoirCouplingMaster&
575AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
576reservoirCouplingMaster_()
577{
578 return *(this->solver_.model().simulator().reservoirCouplingMaster());
579}
580#endif
581
582#ifdef RESERVOIR_COUPLING_ENABLED
583template <class TypeTag>
584template <class Solver>
585ReservoirCouplingSlave&
586AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
587reservoirCouplingSlave_()
588{
589 return *(this->solver_.model().simulator().reservoirCouplingSlave());
590}
591#endif
592
593#ifdef RESERVOIR_COUPLING_ENABLED
594// Description of the reservoir coupling master and slave substep loop
595// -------------------------------------------------------------------
596// The master and slave processes attempts to reach the end of the report step using a series of substeps
597// (also called timesteps). Each substep have an upper limit that is roughly determined by a combination
598// of the keywords TUNING (through the TSINIT, TSMAXZ values), NEXSTEP, WCYCLE, and the start of the
599// next report step (the last substep needs to coincide with this time). Note that NEXTSTEP can be
600// updated from an ACTIONX keyword. Although, this comment focuses on the maximum substep limit,
601// note that there is also a lower limit on the substep length. And the substep sizes will be adjusted
602// automatically (or retried) based on the convergence behavior of the solver and other criteria.
603//
604// When using reservoir coupling, the upper limit on each substep is further limited to: a) not overshoot
605// next report date of a slave reservoir, and b) to keep the flow rate of the slave groups within
606// certain limits. To determine this potential further limitation on the substep, the following procedure
607// is used at the beginning of each master substep:
608// - First, the slaves sends the master the date of their next report step
609// - The master receives the date of the next report step from the slaves
610// - If necessary, the master computes a modified substep length based on the received dates
611// TODO: explain how the substep is limited to keep the flow rate of the slave groups within certain
612// limits. Since this is not implemented yet, this part of the procedure is not explained here.
613//
614// Then, after the master has determined the substep length using the above procedure, it sends it
615// to the slaves. The slaves will now use the end of this substep as a fixed point (like a mini report
616// step), that they will try to reach using a single substep or multiple substeps. The end of the
617// substep received from the master will either conincide with the end of its current report step, or
618// it will end before (it cannot not end after since the master has already adjusted the step to not
619// exceed any report time in a slave). If the end of this substep is earlier in time than its next report
620// date, the slave will enter a loop and wait for the master to send it a new substep.
621// The loop is finished when the end of the received time step conincides with the end of its current
622// report step.
623
624template <class TypeTag>
625template <class Solver>
626SimulatorReport
627AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
628runStepReservoirCouplingMaster_()
629{
630 int iteration = 0;
631 const double original_time_step = this->simulator_timer_.currentStepLength();
632 double current_time{this->simulator_timer_.simulationTimeElapsed()};
633 double step_end_time = current_time + original_time_step;
634 auto current_step_length = original_time_step;
635 SimulatorReport report;
636 while (true) {
637 reservoirCouplingMaster_().receiveNextReportDateFromSlaves();
638 if (iteration == 0) {
639 maybeUpdateTuning_(current_time, current_step_length, /*substep=*/0);
640 }
641 current_step_length = reservoirCouplingMaster_().maybeChopSubStep(
642 current_step_length, current_time);
643 reservoirCouplingMaster_().sendNextTimeStepToSlaves(current_step_length);
644 if (iteration == 0) {
645 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
646 }
647 AdaptiveSimulatorTimer substep_timer{
648 this->simulator_timer_.startDateTime(),
649 /*stepLength=*/current_step_length,
650 /*elapsedTime=*/current_time,
651 /*timeStepEstimate=*/suggestedNextTimestep_(),
652 /*reportStep=*/this->simulator_timer_.reportStepNum(),
653 maxTimeStep_()
654 };
655 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
656 current_time + current_step_length, step_end_time
657 );
658 SubStepIteration<Solver> substepIteration{*this, substep_timer, current_step_length, final_step};
659 const auto sub_steps_report = substepIteration.run();
660 report += sub_steps_report;
661 current_time += current_step_length;
662 if (final_step) {
663 break;
664 }
665 iteration++;
666 }
667 return report;
668}
669#endif
670
671#ifdef RESERVOIR_COUPLING_ENABLED
672template <class TypeTag>
673template <class Solver>
674SimulatorReport
675AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
676runStepReservoirCouplingSlave_()
677{
678 int iteration = 0;
679 const double original_time_step = this->simulator_timer_.currentStepLength();
680 double current_time{this->simulator_timer_.simulationTimeElapsed()};
681 double step_end_time = current_time + original_time_step;
682 SimulatorReport report;
683 while (true) {
684 reservoirCouplingSlave_().sendNextReportDateToMasterProcess();
685 const auto timestep = reservoirCouplingSlave_().receiveNextTimeStepFromMaster();
686 if (iteration == 0) {
687 maybeUpdateTuning_(current_time, original_time_step, /*substep=*/0);
688 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
689 }
690 AdaptiveSimulatorTimer substep_timer{
691 this->simulator_timer_.startDateTime(),
692 /*step_length=*/timestep,
693 /*elapsed_time=*/current_time,
694 /*time_step_estimate=*/suggestedNextTimestep_(),
695 this->simulator_timer_.reportStepNum(),
696 maxTimeStep_()
697 };
698 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
699 current_time + timestep, step_end_time
700 );
701 SubStepIteration<Solver> substepIteration{*this, substep_timer, timestep, final_step};
702 const auto sub_steps_report = substepIteration.run();
703 report += sub_steps_report;
704 current_time += timestep;
705 if (final_step) {
706 break;
707 }
708 iteration++;
709 }
710 return report;
711}
712
713#endif
714
715template <class TypeTag>
716template <class Solver>
717double
718AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
719suggestedNextTimestep_() const
720{
721 return this->adaptive_time_stepping_.suggestedNextStep();
722}
723
724
725
726/************************************************
727 * Private class SubStepIteration public methods
728 ************************************************/
729
730template<class TypeTag>
731template<class Solver>
732AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
733SubStepIteration(SubStepper<Solver>& substepper,
734 AdaptiveSimulatorTimer& substep_timer,
735 const double original_time_step,
736 bool final_step)
737 : substepper_{substepper}
738 , substep_timer_{substep_timer}
739 , original_time_step_{original_time_step}
740 , final_step_{final_step}
741 , adaptive_time_stepping_{substepper.getAdaptiveTimerStepper()}
742{
743}
744
745template <class TypeTag>
746template <class Solver>
747SimulatorReport
748AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
749run()
750{
751 auto& simulator = solver_().model().simulator();
752 auto& problem = simulator.problem();
753 // counter for solver restarts
754 int restarts = 0;
755 SimulatorReport report;
756
757 // sub step time loop
758 while (!this->substep_timer_.done()) {
759 // if we just chopped the timestep due to convergence i.e. restarts>0
760 // we dont want to update the next timestep based on Tuning
761 if (restarts == 0) {
762 maybeUpdateTuningAndTimeStep_();
763 }
764 const double dt = this->substep_timer_.currentStepLength();
765 if (timeStepVerbose_()) {
766 detail::logTimer(this->substep_timer_);
767 }
768
769 const auto substep_report = runSubStep_();
770
771 //Pass substep to eclwriter for summary output
772 problem.setSubStepReport(substep_report);
773 auto& full_report = adaptive_time_stepping_.report();
774 full_report += substep_report;
775 problem.setSimulationReport(full_report);
776
777 report += substep_report;
778
779 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
780 ++this->substep_timer_; // advance by current dt
781
782 const int iterations = getNumIterations_(substep_report);
783 auto dt_estimate = timeStepControlComputeEstimate_(
784 dt, iterations, this->substep_timer_);
785
786 assert(dt_estimate > 0);
787 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
788 restarts = 0; // solver converged, reset restarts counter
789
790 maybeReportSubStep_(substep_report);
791 if (this->final_step_ && this->substep_timer_.done()) {
792 // if the time step is done we do not need to write it as this will be done
793 // by the simulator anyway.
794 }
795 else {
796 report.success.output_write_time += writeOutput_();
797 }
798
799 // set new time step length
800 setTimeStep_(dt_estimate);
801
802 report.success.converged = this->substep_timer_.done();
803 this->substep_timer_.setLastStepFailed(false);
804 }
805 else { // in case of no convergence or time step tolerance test failure
806 this->substep_timer_.setLastStepFailed(true);
807 checkTimeStepMaxRestartLimit_(restarts);
808
809 double new_time_step = restartFactor_() * dt;
810 if (substep_report.time_step_rejected) {
811 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>();
812 const double safetyFactor = Parameters::Get<Parameters::TimeStepControlSafetyFactor>();
813 const double temp_time_step = std::sqrt(safetyFactor * tol / solver_().model().relativeChange()) * dt;
814 if (temp_time_step < dt) { // added in case suggested time step is not a reduction
815 new_time_step = temp_time_step;
816 }
817 }
818 checkTimeStepMinLimit_(new_time_step);
819 bool wells_shut = false;
820 if (new_time_step > minTimeStepBeforeClosingWells_()) {
821 chopTimeStep_(new_time_step);
822 } else {
823 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
824 }
825 if (wells_shut) {
826 setTimeStep_(dt); // retry the old timestep
827 }
828 else {
829 restarts++; // only increase if no wells were shut
830 }
831 }
832 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
833 }
834 updateSuggestedNextStep_();
835 return report;
836}
837
838
839/************************************************
840 * Private class SubStepIteration private methods
841 ************************************************/
842
843
844template<class TypeTag>
845template<class Solver>
846bool
847AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
848checkContinueOnUnconvergedSolution_(double dt) const
849{
850 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
851 if (continue_on_uncoverged_solution && solverVerbose_()) {
852 // NOTE: This method is only called if the solver failed to converge.
853 const auto msg = fmt::format(
854 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
855 "which is the minimum threshold given by option --solver-min-time-step\n",
856 dt, minTimeStep_()
857 );
858 OpmLog::problem(msg);
859 }
860 return continue_on_uncoverged_solution;
861}
862
863template<class TypeTag>
864template<class Solver>
865void
866AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
867checkTimeStepMaxRestartLimit_(const int restarts) const
868{
869 // If we have restarted (i.e. cut the timestep) too
870 // many times, we have failed and throw an exception.
871 if (restarts >= solverRestartMax_()) {
872 const auto msg = fmt::format(
873 "Solver failed to converge after cutting timestep {} times.", restarts
874 );
875 if (solverVerbose_()) {
876 OpmLog::error(msg);
877 }
878 // Use throw directly to prevent file and line
879 throw TimeSteppingBreakdown{msg};
880 }
881}
882
883template<class TypeTag>
884template<class Solver>
885void
886AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
887checkTimeStepMinLimit_(const double new_time_step) const
888{
889 using Meas = UnitSystem::measure;
890 // If we have restarted (i.e. cut the timestep) too
891 // much, we have failed and throw an exception.
892 if (new_time_step < minTimeStep_()) {
893 auto msg = fmt::format("Solver failed to converge after cutting timestep to ");
894 if (Parameters::Get<Parameters::EnableTuning>()) {
895 const UnitSystem& unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
896 msg += fmt::format(
897 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
898 unit_system.from_si(Meas::time, minTimeStep_()),
899 unit_system.name(Meas::time)
900 );
901 }
902 else {
903 msg += fmt::format(
904 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
905 minTimeStep_() / 86400.0
906 );
907 }
908 if (solverVerbose_()) {
909 OpmLog::error(msg);
910 }
911 // Use throw directly to prevent file and line
912 throw TimeSteppingBreakdown{msg};
913 }
914}
915
916template<class TypeTag>
917template<class Solver>
918void
919AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
920chopTimeStep_(const double new_time_step)
921{
922 setTimeStep_(new_time_step);
923 if (solverVerbose_()) {
924 const auto msg = fmt::format("{}\nTimestep chopped to {} days\n",
925 this->cause_of_failure_,
926 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
927 OpmLog::problem(msg);
928 }
929}
930
931template<class TypeTag>
932template<class Solver>
933bool
934AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
935chopTimeStepOrCloseFailingWells_(const double new_time_step)
936{
937 bool wells_shut = false;
938 // We are below the threshold, and will check if there are any
939 // wells that fails repeatedly (that means that it fails in the last three steps)
940 // we should close rather than chopping again.
941 // If we already have chopped the timestep two times that is
942 // new_time_step < minTimeStepBeforeClosingWells_()*restartFactor_()*restartFactor_()
943 // We also shut wells that fails only on this step.
944 const bool requireRepeatedFailures =
945 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
946 const std::set<std::string> failing_wells =
947 detail::consistentlyFailingWells(solver_().model().stepReports(), requireRepeatedFailures);
948
949 if (failing_wells.empty()) {
950 // Found no wells to close, chop the timestep
951 chopTimeStep_(new_time_step);
952 } else {
953 // Close all consistently failing wells that are not under group control
954 std::vector<std::string> shut_wells;
955 for (const auto& well : failing_wells) {
956 const bool was_shut =
957 solver_().model().wellModel().forceShutWellByName(well,
958 this->substep_timer_.simulationTimeElapsed(),
959 /*dont_shut_grup_wells =*/ true);
960 if (was_shut) {
961 shut_wells.push_back(well);
962 }
963 }
964 // If no wells are closed we also try to shut wells under group control
965 if (shut_wells.empty()) {
966 for (const auto& well : failing_wells) {
967 const bool was_shut =
968 solver_().model().wellModel().forceShutWellByName(well,
969 this->substep_timer_.simulationTimeElapsed(),
970 /*dont_shut_grup_wells =*/ false);
971 if (was_shut) {
972 shut_wells.push_back(well);
973 }
974 }
975 }
976 // If still no wells are closed we must fall back to chopping again
977 if (shut_wells.empty()) {
978 chopTimeStep_(new_time_step);
979 } else {
980 wells_shut = true;
981 if (solverVerbose_()) {
982 const std::string msg =
983 fmt::format("\nProblematic well(s) were shut: {}"
984 "(retrying timestep)\n",
985 fmt::join(shut_wells, " "));
986 OpmLog::problem(msg);
987 }
988 }
989 }
990 return wells_shut;
991}
992
993template<class TypeTag>
994template<class Solver>
995boost::posix_time::ptime
996AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
997currentDateTime_() const
998{
999 return simulatorTimer_().currentDateTime();
1000}
1001
1002template<class TypeTag>
1003template<class Solver>
1004int
1005AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1006getNumIterations_(const SimulatorReportSingle &substep_report) const
1007{
1008 if (useNewtonIteration_()) {
1009 return substep_report.total_newton_iterations;
1010 }
1011 else {
1012 return substep_report.total_linear_iterations;
1013 }
1014}
1015
1016template<class TypeTag>
1017template<class Solver>
1018double
1019AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1020growthFactor_() const
1021{
1022 return this->adaptive_time_stepping_.growth_factor_;
1023}
1024
1025template<class TypeTag>
1026template<class Solver>
1027bool
1028AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1029ignoreConvergenceFailure_() const
1030{
1031 return adaptive_time_stepping_.ignore_convergence_failure_;
1032}
1033
1034template<class TypeTag>
1035template<class Solver>
1036double
1037AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1038maxGrowth_() const
1039{
1040 return this->adaptive_time_stepping_.max_growth_;
1041}
1042
1043template<class TypeTag>
1044template<class Solver>
1045void
1046AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1047maybeReportSubStep_(SimulatorReportSingle substep_report) const
1048{
1049 if (timeStepVerbose_()) {
1050 std::ostringstream ss;
1051 substep_report.reportStep(ss);
1052 OpmLog::info(ss.str());
1053 }
1054}
1055
1056template<class TypeTag>
1057template<class Solver>
1058double
1059AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1060maybeRestrictTimeStepGrowth_(const double dt, double dt_estimate, const int restarts) const
1061{
1062 // limit the growth of the timestep size by the growth factor
1063 dt_estimate = std::min(dt_estimate, double(maxGrowth_() * dt));
1064 assert(dt_estimate > 0);
1065 // further restrict time step size growth after convergence problems
1066 if (restarts > 0) {
1067 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1068 }
1069
1070 return dt_estimate;
1071}
1072
1073// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
1074// It has to be called for each substep since TUNING might have been changed for next sub step due
1075// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
1076template<class TypeTag>
1077template<class Solver>
1078void
1079AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1080maybeUpdateTuningAndTimeStep_()
1081{
1082 // TODO: This function is currently only called if NEXTSTEP is activated from ACTIONX or
1083 // if the WCYCLE keyword needs to modify the current timestep. So this method should rather
1084 // be named maybeUpdateTimeStep_() or similar, since it should not update the tuning. However,
1085 // the current definition of the maybeUpdateTuning_() callback is actually calling
1086 // adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning) which is updating the tuning
1087 // see SimulatorFullyImplicitBlackoil::runStep() for more details.
1088 const auto old_value = suggestedNextTimestep_();
1089 if (this->substepper_.maybeUpdateTuning_(this->substep_timer_.simulationTimeElapsed(),
1090 this->substep_timer_.currentStepLength(),
1091 this->substep_timer_.currentStepNum()))
1092 {
1093 // Either NEXTSTEP and WCYCLE wants to change the current time step, but they cannot
1094 // change the current time step directly. Instead, they change the suggested next time step
1095 // by calling updateNEXTSTEP() via the maybeUpdateTuning() callback. We now need to update
1096 // the current time step to the new suggested time step and reset the suggested time step
1097 // to the old value.
1098 setTimeStep_(suggestedNextTimestep_());
1099 setSuggestedNextStep_(old_value);
1100 }
1101}
1102
1103template<class TypeTag>
1104template<class Solver>
1105double
1106AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1107minTimeStepBeforeClosingWells_() const
1108{
1109 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1110}
1111
1112template<class TypeTag>
1113template<class Solver>
1114double
1115AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1116minTimeStep_() const
1117{
1118 return this->adaptive_time_stepping_.min_time_step_;
1119}
1120
1121template<class TypeTag>
1122template<class Solver>
1123double
1124AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1125restartFactor_() const
1126{
1127 return this->adaptive_time_stepping_.restart_factor_;
1128}
1129
1130template<class TypeTag>
1131template<class Solver>
1132SimulatorReportSingle
1133AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1134runSubStep_()
1135{
1136 SimulatorReportSingle substep_report;
1137
1138 auto handleFailure = [this, &substep_report]
1139 (const std::string& failure_reason, const std::exception& e, bool log_exception = true)
1140 {
1141 substep_report = solver_().failureReport();
1142 this->cause_of_failure_ = failure_reason;
1143 if (log_exception && solverVerbose_()) {
1144 OpmLog::debug(std::string("Caught Exception: ") + e.what());
1145 }
1146 };
1147
1148 try {
1149 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1150 if (solverVerbose_()) {
1151 // report number of linear iterations
1152 OpmLog::debug("Overall linear iterations used: "
1153 + std::to_string(substep_report.total_linear_iterations));
1154 }
1155 }
1156 catch (const TooManyIterations& e) {
1157 handleFailure("Solver convergence failure - Iteration limit reached", e);
1158 }
1159 catch (const TimeSteppingBreakdown& e) {
1160 handleFailure(e.what(), e);
1161 }
1162 catch (const ConvergenceMonitorFailure& e) {
1163 handleFailure("Convergence monitor failure", e, /*log_exception=*/false);
1164 }
1165 catch (const LinearSolverProblem& e) {
1166 handleFailure("Linear solver convergence failure", e);
1167 }
1168 catch (const NumericalProblem& e) {
1169 handleFailure("Solver convergence failure - Numerical problem encountered", e);
1170 }
1171 catch (const std::runtime_error& e) {
1172 handleFailure("Runtime error encountered", e);
1173 }
1174 catch (const Dune::ISTLError& e) {
1175 handleFailure("ISTL error - Time step too large", e);
1176 }
1177 catch (const Dune::MatrixBlockError& e) {
1178 handleFailure("Matrix block error", e);
1179 }
1180
1181 return substep_report;
1182}
1183
1184template<class TypeTag>
1185template<class Solver>
1186void
1187AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1188setTimeStep_(double dt_estimate)
1189{
1190 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1191}
1192
1193template<class TypeTag>
1194template<class Solver>
1195Solver&
1196AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1197solver_() const
1198{
1199 return this->substepper_.solver_;
1200}
1201
1202
1203template<class TypeTag>
1204template<class Solver>
1205int
1206AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1207solverRestartMax_() const
1208{
1209 return this->adaptive_time_stepping_.solver_restart_max_;
1210}
1211
1212template<class TypeTag>
1213template<class Solver>
1214void
1215AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1216setSuggestedNextStep_(double step)
1217{
1218 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1219}
1220
1221template <class TypeTag>
1222template <class Solver>
1223const SimulatorTimer&
1224AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1225simulatorTimer_() const
1226{
1227 return this->substepper_.simulator_timer_;
1228}
1229
1230template <class TypeTag>
1231template <class Solver>
1232bool
1233AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1234solverVerbose_() const
1235{
1236 return this->adaptive_time_stepping_.solver_verbose_;
1237}
1238
1239template<class TypeTag>
1240template<class Solver>
1241boost::posix_time::ptime
1242AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1243startDateTime_() const
1244{
1245 return simulatorTimer_().startDateTime();
1246}
1247
1248template <class TypeTag>
1249template <class Solver>
1250double
1251AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1252suggestedNextTimestep_() const
1253{
1254 return this->adaptive_time_stepping_.suggestedNextStep();
1255}
1256
1257template <class TypeTag>
1258template <class Solver>
1259double
1260AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1261timeStepControlComputeEstimate_(const double dt, const int iterations,
1262 const AdaptiveSimulatorTimer& substepTimer) const
1263{
1264 // create object to compute the time error, simply forwards the call to the model
1265 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1266 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1267 dt, iterations, relative_change, substepTimer);
1268}
1269
1270template <class TypeTag>
1271template <class Solver>
1272bool
1273AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1274timeStepVerbose_() const
1275{
1276 return this->adaptive_time_stepping_.timestep_verbose_;
1277}
1278
1279// The suggested time step is the stepsize that will be used as a first try for
1280// the next sub step. It is updated at the end of each substep. It can also be
1281// updated by the TUNING or NEXTSTEP keywords at the beginning of each report step or
1282// at the beginning of each substep by the ACTIONX keyword (via NEXTSTEP), this is
1283// done by the maybeUpdateTuning_() method which is called at the beginning of each substep
1284// (and the begginning of each report step). Note that the WCYCLE keyword can also update the
1285// suggested time step via the maybeUpdateTuning_() method.
1286template <class TypeTag>
1287template <class Solver>
1288void
1289AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1290updateSuggestedNextStep_()
1291{
1292 auto suggested_next_step = this->substep_timer_.currentStepLength();
1293 if (! std::isfinite(suggested_next_step)) { // check for NaN
1294 suggested_next_step = this->original_time_step_;
1295 }
1296 if (timeStepVerbose_()) {
1297 std::ostringstream ss;
1298 this->substep_timer_.report(ss);
1299 ss << "Suggested next step size = "
1300 << unit::convert::to(suggested_next_step, unit::day) << " (days)" << std::endl;
1301 OpmLog::debug(ss.str());
1302 }
1303 setSuggestedNextStep_(suggested_next_step);
1304}
1305
1306template <class TypeTag>
1307template <class Solver>
1308bool
1309AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1310useNewtonIteration_() const
1311{
1312 return this->adaptive_time_stepping_.use_newton_iteration_;
1313}
1314
1315template <class TypeTag>
1316template <class Solver>
1317double
1318AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1319writeOutput_() const
1320{
1321 time::StopWatch perf_timer;
1322 perf_timer.start();
1323 auto& problem = solver_().model().simulator().problem();
1324 problem.writeOutput(true);
1325 return perf_timer.secsSinceStart();
1326}
1327
1328/************************************************
1329 * Private class SolutionTimeErrorSolverWrapper
1330 * **********************************************/
1331
1332template<class TypeTag>
1333template<class Solver>
1334AdaptiveTimeStepping<TypeTag>::
1335SolutionTimeErrorSolverWrapper<Solver>::
1336SolutionTimeErrorSolverWrapper(const Solver& solver)
1337 : solver_{solver}
1338{}
1339
1340template<class TypeTag>
1341template<class Solver>
1342double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange() const
1343{
1344 // returns: || u^n+1 - u^n || / || u^n+1 ||
1345 return solver_.model().relativeChange();
1346}
1347
1348} // namespace Opm
1349
1350#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:427
void setSuggestedNextStep(const double x)
Definition: AdaptiveTimeStepping_impl.hpp:298
double suggestedNextStep() const
Definition: AdaptiveTimeStepping_impl.hpp:306
std::function< bool(const double, const double, const int)> TuningUpdateCallback
Definition: AdaptiveTimeStepping.hpp:80
static AdaptiveTimeStepping< TypeTag > serializationTestObjectSimple()
Definition: AdaptiveTimeStepping_impl.hpp:281
bool ignore_convergence_failure_
continue instead of stop when minimum time step is reached
Definition: AdaptiveTimeStepping.hpp:255
void serializeOp(Serializer &serializer)
Definition: AdaptiveTimeStepping_impl.hpp:210
void updateTUNING(double max_next_tstep, const Tuning &tuning)
Definition: AdaptiveTimeStepping_impl.hpp:335
TimeStepControlType time_step_control_type_
type of time step control object
Definition: AdaptiveTimeStepping.hpp:248
const TimeStepControlInterface & timeStepControl() const
Definition: AdaptiveTimeStepping_impl.hpp:314
bool full_timestep_initially_
beginning with the size of the time step from data file
Definition: AdaptiveTimeStepping.hpp:260
bool operator==(const AdaptiveTimeStepping< TypeTag > &rhs)
Definition: AdaptiveTimeStepping_impl.hpp:134
SimulatorReport step(const SimulatorTimer &simulator_timer, Solver &solver, const bool is_event, const TuningUpdateCallback &tuning_updater)
step method that acts like the solver::step method in a sub cycle of time steps
Definition: AdaptiveTimeStepping_impl.hpp:195
double growth_factor_
factor to multiply time step when solver recovered from failed convergence
Definition: AdaptiveTimeStepping.hpp:251
double restart_factor_
factor to multiply time step with when solver fails to converge
Definition: AdaptiveTimeStepping.hpp:250
double min_time_step_
minimal allowed time step size before throwing
Definition: AdaptiveTimeStepping.hpp:254
void updateNEXTSTEP(double max_next_tstep)
Definition: AdaptiveTimeStepping_impl.hpp:323
static AdaptiveTimeStepping< TypeTag > serializationTestObjectHardcoded()
Definition: AdaptiveTimeStepping_impl.hpp:257
TimeStepController time_step_control_
time step control object
Definition: AdaptiveTimeStepping.hpp:249
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPIDIt()
Definition: AdaptiveTimeStepping_impl.hpp:273
double min_time_step_before_shutting_problematic_wells_
< shut problematic wells when time step size in days are less than this
Definition: AdaptiveTimeStepping.hpp:265
static AdaptiveTimeStepping< TypeTag > serializationTestObject3rdOrder()
Definition: AdaptiveTimeStepping_impl.hpp:289
SimulatorReport & report()
Definition: AdaptiveTimeStepping_impl.hpp:249
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPID()
Definition: AdaptiveTimeStepping_impl.hpp:265
static void registerParameters()
Definition: AdaptiveTimeStepping_impl.hpp:180
bool use_newton_iteration_
use newton iteration count for adaptive time step control
Definition: AdaptiveTimeStepping.hpp:262
Definition: SimulatorTimer.hpp:39
Definition: TimeStepControlInterface.hpp:51
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:185
void logTimer(const AdaptiveSimulatorTimer &substep_timer)
void registerAdaptiveParameters()
std::set< std::string > consistentlyFailingWells(const std::vector< StepReport > &sr, bool requireRepeatedFailures)
std::tuple< TimeStepControlType, std::unique_ptr< TimeStepControlInterface >, bool > createController(const UnitSystem &unitSystem)
Definition: blackoilboundaryratevector.hh:39
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
This file provides the infrastructure to retrieve run-time parameters.
static bool compare_gt_or_eq(double a, double b)
Determines if a is greater than b within the specified tolerance.
Definition: SimulatorReport.hpp:122