opm-simulators
Loading...
Searching...
No Matches
AdaptiveTimeStepping_impl.hpp
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>
26#include <opm/simulators/timestepping/AdaptiveTimeStepping.hpp>
27#include <opm/simulators/timestepping/AdaptiveSimulatorTimer.hpp>
28#endif
29
30#include <dune/common/timer.hh>
31#include <dune/istl/istlexception.hh>
32
33#include <opm/common/Exceptions.hpp>
34#include <opm/common/ErrorMacros.hpp>
35#include <opm/common/OpmLog/OpmLog.hpp>
36
37#include <opm/grid/utility/StopWatch.hpp>
38
39#include <opm/input/eclipse/Schedule/Tuning.hpp>
40
41#include <opm/input/eclipse/Units/Units.hpp>
42#include <opm/input/eclipse/Units/UnitSystem.hpp>
43
45
46#include <opm/simulators/timestepping/EclTimeSteppingParams.hpp>
47
48#include <algorithm>
49#include <cassert>
50#include <cmath>
51#include <sstream>
52#include <stdexcept>
53
54#include <boost/date_time/posix_time/posix_time.hpp>
55#include <fmt/format.h>
56#include <fmt/ranges.h>
57namespace Opm {
58/*********************************************
59 * Public methods of AdaptiveTimeStepping
60 * ******************************************/
61
62
64template<class TypeTag>
65AdaptiveTimeStepping<TypeTag>::
66AdaptiveTimeStepping(const UnitSystem& unit_system,
67 const SimulatorReport& report,
68 const double max_next_tstep,
69 const bool terminal_output
70)
72 , restart_factor_{Parameters::Get<Parameters::SolverRestartFactor<Scalar>>()} // 0.33
73 , growth_factor_{Parameters::Get<Parameters::SolverGrowthFactor<Scalar>>()} // 2.0
74 , max_growth_{Parameters::Get<Parameters::SolverMaxGrowth<Scalar>>()} // 3.0
76 Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60} // 365.25
78 unit_system.to_si(UnitSystem::measure::time,
79 Parameters::Get<Parameters::SolverMinTimeStep<Scalar>>())} // 1e-12;
81 Parameters::Get<Parameters::SolverContinueOnConvergenceFailure>()} // false;
82 , solver_restart_max_{Parameters::Get<Parameters::SolverMaxRestarts>()} // 10
83 , solver_verbose_{Parameters::Get<Parameters::SolverVerbosity>() > 0 && terminal_output} // 2
84 , timestep_verbose_{Parameters::Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output} // 2
86 (max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>()
87 : max_next_tstep) * 24 * 60 * 60} // 1.0
88 , full_timestep_initially_{Parameters::Get<Parameters::FullTimeStepInitially>()} // false
90 Parameters::Get<Parameters::TimeStepAfterEventInDays<Scalar>>() * 24 * 60 * 60} // 1e30
93 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
94 , report_(report)
95{
96 init_(unit_system);
97}
98
105template<class TypeTag>
106AdaptiveTimeStepping<TypeTag>::
107AdaptiveTimeStepping(double max_next_tstep,
108 const Tuning& tuning,
109 const UnitSystem& unit_system,
110 const SimulatorReport& report,
111 const bool terminal_output
112)
114 , restart_factor_{tuning.TSFCNV}
115 , growth_factor_{tuning.TFDIFF}
116 , max_growth_{tuning.TSFMAX}
117 , max_time_step_{tuning.TSMAXZ} // 365.25
118 , min_time_step_{tuning.TSMINZ} // 0.1;
120 , solver_restart_max_{Parameters::Get<Parameters::SolverMaxRestarts>()} // 10
121 , solver_verbose_{Parameters::Get<Parameters::SolverVerbosity>() > 0 && terminal_output} // 2
122 , timestep_verbose_{Parameters::Get<Parameters::TimeStepVerbosity>() > 0 && terminal_output} // 2
124 max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>() * 24 * 60 * 60
125 : max_next_tstep} // 1.0
126 , full_timestep_initially_{Parameters::Get<Parameters::FullTimeStepInitially>()} // false
127 , timestep_after_event_{tuning.TMAXWC} // 1e30
128 , use_newton_iteration_{false}
130 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
131 , report_(report)
132{
133 init_(unit_system);
134}
135
136template<class TypeTag>
137bool
140{
141 if (this->time_step_control_type_ != rhs.time_step_control_type_ ||
142 (this->time_step_control_ && !rhs.time_step_control_) ||
143 (!this->time_step_control_ && rhs.time_step_control_)) {
144 return false;
145 }
146
147 bool result = false;
148 switch (this->time_step_control_type_) {
149 case TimeStepControlType::HardCodedTimeStep:
150 result = castAndComp<HardcodedTimeStepControl>(rhs);
151 break;
152 case TimeStepControlType::PIDAndIterationCount:
153 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
154 break;
155 case TimeStepControlType::SimpleIterationCount:
156 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
157 break;
158 case TimeStepControlType::PID:
159 result = castAndComp<PIDTimeStepControl>(rhs);
160 break;
161 case TimeStepControlType::General3rdOrder:
162 result = castAndComp<General3rdOrderController>(rhs);
163 break;
164 }
165
166 return result &&
167 this->restart_factor_ == rhs.restart_factor_ &&
168 this->growth_factor_ == rhs.growth_factor_ &&
169 this->max_growth_ == rhs.max_growth_ &&
170 this->max_time_step_ == rhs.max_time_step_ &&
171 this->min_time_step_ == rhs.min_time_step_ &&
172 this->ignore_convergence_failure_ == rhs.ignore_convergence_failure_ &&
173 this->solver_restart_max_== rhs.solver_restart_max_ &&
174 this->solver_verbose_ == rhs.solver_verbose_ &&
175 this->full_timestep_initially_ == rhs.full_timestep_initially_ &&
176 this->timestep_after_event_ == rhs.timestep_after_event_ &&
177 this->use_newton_iteration_ == rhs.use_newton_iteration_ &&
178 this->min_time_step_before_shutting_problematic_wells_ ==
180}
181
182template<class TypeTag>
183void
186{
187 registerEclTimeSteppingParameters<Scalar>();
188 detail::registerAdaptiveParameters();
189}
190
199template<class TypeTag>
200template <class Solver>
203step(const SimulatorTimer& simulator_timer,
204 Solver& solver,
205 const bool is_event,
206 const TuningUpdateCallback& tuning_updater)
207{
208 SubStepper<Solver> sub_stepper{
209 *this, simulator_timer, solver, is_event, tuning_updater,
210 };
211 return sub_stepper.run();
212}
213
214template<class TypeTag>
215template<class Serializer>
216void
218serializeOp(Serializer& serializer)
219{
220 serializer(this->time_step_control_type_);
221 switch (this->time_step_control_type_) {
222 case TimeStepControlType::HardCodedTimeStep:
223 allocAndSerialize<HardcodedTimeStepControl>(serializer);
224 break;
225 case TimeStepControlType::PIDAndIterationCount:
226 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
227 break;
228 case TimeStepControlType::SimpleIterationCount:
229 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
230 break;
231 case TimeStepControlType::PID:
232 allocAndSerialize<PIDTimeStepControl>(serializer);
233 break;
234 case TimeStepControlType::General3rdOrder:
235 allocAndSerialize<General3rdOrderController>(serializer);
236 break;
237 }
238 serializer(this->restart_factor_);
239 serializer(this->growth_factor_);
240 serializer(this->max_growth_);
241 serializer(this->max_time_step_);
242 serializer(this->min_time_step_);
243 serializer(this->ignore_convergence_failure_);
244 serializer(this->solver_restart_max_);
245 serializer(this->solver_verbose_);
246 serializer(this->timestep_verbose_);
247 serializer(this->suggested_next_timestep_);
248 serializer(this->full_timestep_initially_);
249 serializer(this->timestep_after_event_);
250 serializer(this->use_newton_iteration_);
251 serializer(this->min_time_step_before_shutting_problematic_wells_);
252}
253
254template<class TypeTag>
257report()
258{
259 return report_;
260}
261
262template<class TypeTag>
266{
267 return serializationTestObject_<HardcodedTimeStepControl>();
268}
269
270template<class TypeTag>
274{
275 return serializationTestObject_<PIDTimeStepControl>();
276}
277
278template<class TypeTag>
282{
283 return serializationTestObject_<PIDAndIterationCountTimeStepControl>();
284}
285
286template<class TypeTag>
290{
291 return serializationTestObject_<SimpleIterationCountTimeStepControl>();
292}
293
294template<class TypeTag>
298{
299 return serializationTestObject_<General3rdOrderController>();
300}
301
302
303template<class TypeTag>
304void
306setSuggestedNextStep(const double x)
307{
308 this->suggested_next_timestep_ = x;
309}
310
311template<class TypeTag>
312double
314suggestedNextStep() const
315{
316 return this->suggested_next_timestep_;
317}
318
319template<class TypeTag>
322timeStepControl() const
323{
324 return *this->time_step_control_;
325}
326
327
328template<class TypeTag>
329void
331updateNEXTSTEP(double max_next_tstep)
332{
333 // \Note Only update next suggested step if TSINIT was explicitly
334 // set in TUNING or NEXTSTEP is active.
335 if (max_next_tstep > 0) {
336 this->suggested_next_timestep_ = max_next_tstep;
337 }
338}
339
340template<class TypeTag>
341void
343updateTUNING(double max_next_tstep, const Tuning& tuning)
344{
345 this->restart_factor_ = tuning.TSFCNV;
346 this->growth_factor_ = tuning.TFDIFF;
347 this->max_growth_ = tuning.TSFMAX;
348 this->max_time_step_ = tuning.TSMAXZ;
349 updateNEXTSTEP(max_next_tstep);
350 this->timestep_after_event_ = tuning.TMAXWC;
351}
352
353/*********************************************
354 * Private methods of AdaptiveTimeStepping
355 * ******************************************/
356
357template<class TypeTag>
358template<class T, class Serializer>
359void
361allocAndSerialize(Serializer& serializer)
362{
363 if (!serializer.isSerializing()) {
364 this->time_step_control_ = std::make_unique<T>();
365 }
366 serializer(*static_cast<T*>(this->time_step_control_.get()));
367}
368
369template<class TypeTag>
370template<class T>
371bool
374{
375 const T* lhs = static_cast<const T*>(this->time_step_control_.get());
376 const T* rhs = static_cast<const T*>(Rhs.time_step_control_.get());
377 return *lhs == *rhs;
378}
379
380template<class TypeTag>
381void
383maybeModifySuggestedTimeStepAtBeginningOfReportStep_(const double original_time_step,
384 bool is_event)
385{
386 // init last time step as a fraction of the given time step
387 if (this->suggested_next_timestep_ < 0) {
388 this->suggested_next_timestep_ = this->restart_factor_ * original_time_step;
389 }
390
391 if (this->full_timestep_initially_) {
392 this->suggested_next_timestep_ = original_time_step;
393 }
394
395 // use seperate time step after event
396 if (is_event && this->timestep_after_event_ > 0) {
397 this->suggested_next_timestep_ = this->timestep_after_event_;
398 }
399}
400
401template<class TypeTag>
402template<class Controller>
406{
408
409 result.restart_factor_ = 1.0;
410 result.growth_factor_ = 2.0;
411 result.max_growth_ = 3.0;
412 result.max_time_step_ = 4.0;
413 result.min_time_step_ = 5.0;
414 result.ignore_convergence_failure_ = true;
415 result.solver_restart_max_ = 6;
416 result.solver_verbose_ = true;
417 result.timestep_verbose_ = true;
418 result.suggested_next_timestep_ = 7.0;
419 result.full_timestep_initially_ = true;
420 result.use_newton_iteration_ = true;
421 result.min_time_step_before_shutting_problematic_wells_ = 9.0;
422 result.time_step_control_type_ = Controller::Type;
423 result.time_step_control_ =
424 std::make_unique<Controller>(Controller::serializationTestObject());
425
426 return result;
427}
428
429/*********************************************
430 * Protected methods of AdaptiveTimeStepping
431 * ******************************************/
432
433template<class TypeTag>
435init_(const UnitSystem& unitSystem)
436{
437 std::tie(time_step_control_type_,
438 time_step_control_,
439 use_newton_iteration_) = detail::createController(unitSystem);
440 // make sure growth factor is something reasonable
441 if (this->growth_factor_ < 1.0) {
442 OPM_THROW(std::runtime_error,
443 "Growth factor cannot be less than 1.");
444 }
445}
446
447
448
449/************************************************
450 * Private class SubStepper public methods
451 ************************************************/
452
453template<class TypeTag>
454template<class Solver>
456SubStepper(AdaptiveTimeStepping<TypeTag>& adaptive_time_stepping,
457 const SimulatorTimer& simulator_timer,
458 Solver& solver,
459 const bool is_event,
460 const TuningUpdateCallback& tuning_updater)
461 : adaptive_time_stepping_{adaptive_time_stepping}
462 , simulator_timer_{simulator_timer}
463 , solver_{solver}
464 , is_event_{is_event}
465 , tuning_updater_{tuning_updater}
466{
467}
468
469template<class TypeTag>
470template<class Solver>
474{
475 return adaptive_time_stepping_;
476}
477
478template<class TypeTag>
479template<class Solver>
482run()
483{
484#ifdef RESERVOIR_COUPLING_ENABLED
485 if (isReservoirCouplingSlave_() && reservoirCouplingSlave_().activated()) {
486 return runStepReservoirCouplingSlave_();
487 }
488 else if (isReservoirCouplingMaster_() && reservoirCouplingMaster_().activated()) {
489 return runStepReservoirCouplingMaster_();
490 }
491 else {
492 return runStepOriginal_();
493 }
494#else
495 return runStepOriginal_();
496#endif
497}
498
499/************************************************
500 * Private class SubStepper private methods
501 ************************************************/
502
503#ifdef RESERVOIR_COUPLING_ENABLED
504template<class TypeTag>
505template<class Solver>
506bool
509{
510 return this->solver_.model().simulator().reservoirCouplingMaster() != nullptr;
511}
512
513template<class TypeTag>
514template<class Solver>
515bool
518{
519 return this->solver_.model().simulator().reservoirCouplingSlave() != nullptr;
520}
521#endif
522
523template<class TypeTag>
524template<class Solver>
525void
527maybeModifySuggestedTimeStepAtBeginningOfReportStep_(const double original_time_step)
528{
529 this->adaptive_time_stepping_.maybeModifySuggestedTimeStepAtBeginningOfReportStep_(
530 original_time_step, this->is_event_
531 );
532}
533
534// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
535// It has to be called for each substep since TUNING might have been changed for next sub step due
536// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
537template<class TypeTag>
538template<class Solver>
539bool
541maybeUpdateTuning_(double elapsed, double dt, int sub_step_number) const
542{
543 return this->tuning_updater_(elapsed, dt, sub_step_number);
544}
545
546template<class TypeTag>
547template<class Solver>
548double
550maxTimeStep_() const
551{
552 return this->adaptive_time_stepping_.max_time_step_;
553}
554
555template <class TypeTag>
556template <class Solver>
560{
561 const auto elapsed = this->simulator_timer_.simulationTimeElapsed();
562 const auto original_time_step = this->simulator_timer_.currentStepLength();
563 const auto report_step = this->simulator_timer_.reportStepNum();
564 maybeUpdateTuning_(elapsed, original_time_step, report_step);
565 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(original_time_step);
566
567 AdaptiveSimulatorTimer substep_timer{
568 this->simulator_timer_.startDateTime(),
569 original_time_step,
570 elapsed,
571 suggestedNextTimestep_(),
572 report_step,
573 maxTimeStep_()
574 };
575 SubStepIteration<Solver> substepIteration{*this, substep_timer, original_time_step, /*final_step=*/true};
576 return substepIteration.run();
577}
578
579#ifdef RESERVOIR_COUPLING_ENABLED
580template <class TypeTag>
581template <class Solver>
585{
586 return *(this->solver_.model().simulator().reservoirCouplingMaster());
587}
588#endif
589
590#ifdef RESERVOIR_COUPLING_ENABLED
591template <class TypeTag>
592template <class Solver>
596{
597 return *(this->solver_.model().simulator().reservoirCouplingSlave());
598}
599#endif
600
601#ifdef RESERVOIR_COUPLING_ENABLED
602// Description of the reservoir coupling master and slave substep loop
603// -------------------------------------------------------------------
604// The master and slave processes attempts to reach the end of the report step using a series of substeps
605// (also called timesteps). Each substep have an upper limit that is roughly determined by a combination
606// of the keywords TUNING (through the TSINIT, TSMAXZ values), NEXSTEP, WCYCLE, and the start of the
607// next report step (the last substep needs to coincide with this time). Note that NEXTSTEP can be
608// updated from an ACTIONX keyword. Although, this comment focuses on the maximum substep limit,
609// note that there is also a lower limit on the substep length. And the substep sizes will be adjusted
610// automatically (or retried) based on the convergence behavior of the solver and other criteria.
611//
612// When using reservoir coupling, the upper limit on each substep is further limited to: a) not overshoot
613// next report date of a slave reservoir, and b) to keep the flow rate of the slave groups within
614// certain limits. To determine this potential further limitation on the substep, the following procedure
615// is used at the beginning of each master substep:
616// - First, the slaves sends the master the date of their next report step
617// - The master receives the date of the next report step from the slaves
618// - If necessary, the master computes a modified substep length based on the received dates
619// TODO: explain how the substep is limited to keep the flow rate of the slave groups within certain
620// limits. Since this is not implemented yet, this part of the procedure is not explained here.
621//
622// Then, after the master has determined the substep length using the above procedure, it sends it
623// to the slaves. The slaves will now use the end of this substep as a fixed point (like a mini report
624// step), that they will try to reach using a single substep or multiple substeps. The end of the
625// substep received from the master will either conincide with the end of its current report step, or
626// it will end before (it cannot not end after since the master has already adjusted the step to not
627// exceed any report time in a slave). If the end of this substep is earlier in time than its next report
628// date, the slave will enter a loop and wait for the master to send it a new substep.
629// The loop is finished when the end of the received time step conincides with the end of its current
630// report step.
631
632template <class TypeTag>
633template <class Solver>
637{
638 int iteration = 0;
639 const double original_time_step = this->simulator_timer_.currentStepLength();
640 double current_time{this->simulator_timer_.simulationTimeElapsed()};
641 double step_end_time = current_time + original_time_step;
642 auto current_step_length = original_time_step;
643 auto report_step_idx = this->simulator_timer_.currentStepNum();
644 if (report_step_idx == 0 && iteration == 0) {
645 reservoirCouplingMaster_().initTimeStepping();
646 }
647 SimulatorReport report;
648 // The master needs to know which slaves have activated before it can start the substep loop
649 reservoirCouplingMaster_().maybeReceiveActivationHandshakeFromSlaves(current_time);
650 while (true) {
651 reservoirCouplingMaster_().sendDontTerminateSignalToSlaves(); // Tell the slaves to keep running.
652 reservoirCouplingMaster_().receiveNextReportDateFromSlaves();
653 bool start_of_report_step = (iteration == 0);
654 if (start_of_report_step) {
655 reservoirCouplingMaster_().initStartOfReportStep(report_step_idx);
656 }
657 current_step_length = reservoirCouplingMaster_().maybeChopSubStep(
658 current_step_length, current_time);
659 auto num_active = reservoirCouplingMaster_().numActivatedSlaves();
660 OpmLog::info(fmt::format(
661 "\nChoosing next sync time between master and {} active slave {}: {:.2f} days",
662 num_active, (num_active == 1 ? "process" : "processes"),
663 current_step_length / unit::day
664 ));
665 reservoirCouplingMaster_().sendNextTimeStepToSlaves(current_step_length);
666 if (start_of_report_step) {
667 maybeUpdateTuning_(current_time, current_step_length, /*substep=*/0);
668 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
669 }
670 AdaptiveSimulatorTimer substep_timer{
671 this->simulator_timer_.startDateTime(),
672 /*stepLength=*/current_step_length,
673 /*elapsedTime=*/current_time,
674 /*timeStepEstimate=*/suggestedNextTimestep_(),
675 /*reportStep=*/this->simulator_timer_.reportStepNum(),
676 maxTimeStep_()
677 };
678 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
679 current_time + current_step_length, step_end_time
680 );
681 // Mark this as the first substep of the "sync" timestep. This flag controls
682 // whether master-slave data exchange should occur in beginTimeStep() in the well model.
683 // It will be cleared after the first runSubStep_() call.
684 reservoirCouplingMaster_().setFirstSubstepOfSyncTimestep(true);
685 // After the first master substep completes, timeStepSucceeded() will
686 // block until slaves finish the sync step and send production data.
687 // This ensures correct summary output for all subsequent substeps.
688 reservoirCouplingMaster_().setNeedsSlaveDataReceive(true);
689 SubStepIteration<Solver> substepIteration{*this, substep_timer, current_step_length, final_step};
690 const auto sub_steps_report = substepIteration.run();
691 report += sub_steps_report;
692 current_time += current_step_length;
693 if (final_step) {
694 break;
695 }
696 iteration++;
697 }
698 return report;
699}
700#endif
701
702#ifdef RESERVOIR_COUPLING_ENABLED
703template <class TypeTag>
704template <class Solver>
708{
709 int iteration = 0;
710 const double original_time_step = this->simulator_timer_.currentStepLength();
711 double current_time{this->simulator_timer_.simulationTimeElapsed()};
712 double step_end_time = current_time + original_time_step;
713 SimulatorReport report;
714 auto report_step_idx = this->simulator_timer_.currentStepNum();
715 if (report_step_idx == 0 && iteration == 0) {
716 reservoirCouplingSlave_().initTimeStepping();
717 }
718 while (true) {
719 bool start_of_report_step = (iteration == 0);
720 if (reservoirCouplingSlave_().maybeReceiveTerminateSignalFromMaster()) {
721 // Call MPI_Comm_disconnect() to terminate the MPI communicator, etc..
722 break;
723 }
724 reservoirCouplingSlave_().sendNextReportDateToMasterProcess();
725 const auto timestep = reservoirCouplingSlave_().receiveNextTimeStepFromMaster();
726 if (start_of_report_step) {
727 maybeUpdateTuning_(current_time, original_time_step, /*substep=*/0);
728 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
729 }
730 AdaptiveSimulatorTimer substep_timer{
731 this->simulator_timer_.startDateTime(),
732 /*step_length=*/timestep,
733 /*elapsed_time=*/current_time,
734 /*time_step_estimate=*/suggestedNextTimestep_(),
735 this->simulator_timer_.reportStepNum(),
736 maxTimeStep_()
737 };
738 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
739 current_time + timestep, step_end_time
740 );
741 // Mark this as the first substep of the "sync" timestep. This flag controls
742 // whether master-slave data exchange should occur in beginTimeStep() in the well model.
743 // It will be cleared after the first runSubStep_() call.
744 reservoirCouplingSlave_().setFirstSubstepOfSyncTimestep(true);
745 SubStepIteration<Solver> substepIteration{*this, substep_timer, timestep, final_step};
746 const auto sub_steps_report = substepIteration.run();
747 report += sub_steps_report;
748 current_time += timestep;
749 if (final_step) {
750 break;
751 }
752 iteration++;
753 }
754 return report;
755}
756
757#endif
758
759template <class TypeTag>
760template <class Solver>
761double
764{
765 return this->adaptive_time_stepping_.suggestedNextStep();
766}
767
768
769
770/************************************************
771 * Private class SubStepIteration public methods
772 ************************************************/
773
774template<class TypeTag>
775template<class Solver>
777SubStepIteration(SubStepper<Solver>& substepper,
778 AdaptiveSimulatorTimer& substep_timer,
779 const double original_time_step,
780 bool final_step)
781 : substepper_{substepper}
782 , substep_timer_{substep_timer}
783 , original_time_step_{original_time_step}
784 , final_step_{final_step}
785 , adaptive_time_stepping_{substepper.getAdaptiveTimerStepper()}
786{
787}
788
789template <class TypeTag>
790template <class Solver>
793run()
794{
795 auto& simulator = solver_().model().simulator();
796 auto& problem = simulator.problem();
797 // counter for solver restarts
798 int restarts = 0;
799 SimulatorReport report;
800
801 // sub step time loop
802 while (!this->substep_timer_.done()) {
803 // if we just chopped the timestep due to convergence i.e. restarts>0
804 // we dont want to update the next timestep based on Tuning
805 if (restarts == 0) {
806 maybeUpdateTuningAndTimeStep_();
807 }
808 const double dt = this->substep_timer_.currentStepLength();
809 if (timeStepVerbose_()) {
810 detail::logTimer(this->substep_timer_);
811 }
812
813 maybeUpdateLastSubstepOfSyncTimestep_(dt); // Needed for reservoir coupling
814 auto substep_report = runSubStep_();
815 markFirstSubStepAsFinished_(); // Needed for reservoir coupling
816
817 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
818 Dune::Timer perfTimer;
819 perfTimer.start();
820 // Pass substep to eclwriter for summary output
821 problem.setSubStepReport(substep_report);
822 auto& full_report = adaptive_time_stepping_.report();
823 full_report += substep_report;
824 problem.setSimulationReport(full_report);
825 problem.endTimeStep();
826 substep_report.pre_post_time += perfTimer.stop();
827
828 report += substep_report;
829
830 ++this->substep_timer_; // advance by current dt
831
832 const int iterations = getNumIterations_(substep_report);
833 auto dt_estimate = timeStepControlComputeEstimate_(
834 dt, iterations, this->substep_timer_);
835
836 assert(dt_estimate > 0);
837 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
838 restarts = 0; // solver converged, reset restarts counter
839
840 maybeReportSubStep_(substep_report);
841 if (this->final_step_ && this->substep_timer_.done()) {
842 // if the time step is done we do not need to write it as this will be done
843 // by the simulator anyway.
844 }
845 else {
846 report.success.output_write_time += writeOutput_();
847 }
848
849 // set new time step length
850 setTimeStep_(dt_estimate);
851
852 report.success.converged = this->substep_timer_.done();
853 this->substep_timer_.setLastStepFailed(false);
854 }
855 else { // in case of no convergence or time step tolerance test failure
856 report += substep_report;
857 this->substep_timer_.setLastStepFailed(true);
858 checkTimeStepMaxRestartLimit_(restarts);
859
860 double new_time_step = restartFactor_() * dt;
861 if (substep_report.time_step_rejected) {
864 const double temp_time_step = std::sqrt(safetyFactor * tol / solver_().model().relativeChange()) * dt;
865 if (temp_time_step < dt) { // added in case suggested time step is not a reduction
866 new_time_step = temp_time_step;
867 }
868 }
869 checkTimeStepMinLimit_(new_time_step);
870 bool wells_shut = false;
871 if (new_time_step > minTimeStepBeforeClosingWells_()) {
872 chopTimeStep_(new_time_step);
873 } else {
874 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
875 }
876 if (wells_shut) {
877 setTimeStep_(dt); // retry the old timestep
878 }
879 else {
880 restarts++; // only increase if no wells were shut
881 }
882 }
883 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
884 }
885 updateSuggestedNextStep_();
886 return report;
887}
888
889
890/************************************************
891 * Private class SubStepIteration private methods
892 ************************************************/
893
894
895template<class TypeTag>
896template<class Solver>
897bool
900{
901 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
902 if (continue_on_uncoverged_solution && solverVerbose_()) {
903 // NOTE: This method is only called if the solver failed to converge.
904 const auto msg = fmt::format(
905 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
906 "which is the minimum threshold given by option --solver-min-time-step\n",
907 dt, minTimeStep_()
908 );
909 OpmLog::problem(msg);
910 }
911 return continue_on_uncoverged_solution;
912}
913
914template<class TypeTag>
915template<class Solver>
916void
918checkTimeStepMaxRestartLimit_(const int restarts) const
919{
920 // If we have restarted (i.e. cut the timestep) too
921 // many times, we have failed and throw an exception.
922 if (restarts >= solverRestartMax_()) {
923 const auto msg = fmt::format(
924 fmt::runtime("Solver failed to converge after cutting timestep {} times."), restarts
925 );
926 if (solverVerbose_()) {
927 OpmLog::error(msg);
928 }
929 // Use throw directly to prevent file and line
930 throw TimeSteppingBreakdown{msg};
931 }
932}
933
934template<class TypeTag>
935template<class Solver>
936void
938checkTimeStepMinLimit_(const double new_time_step) const
939{
940 using Meas = UnitSystem::measure;
941 // If we have restarted (i.e. cut the timestep) too
942 // much, we have failed and throw an exception.
943 if (new_time_step < minTimeStep_()) {
944 std::string msg = "Solver failed to converge after cutting timestep to ";
946 const UnitSystem& unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
947 msg += fmt::format(
948 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
949 unit_system.from_si(Meas::time, minTimeStep_()),
950 unit_system.name(Meas::time)
951 );
952 }
953 else {
954 msg += fmt::format(
955 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
956 minTimeStep_() / 86400.0
957 );
958 }
959 if (solverVerbose_()) {
960 OpmLog::error(msg);
961 }
962 // Use throw directly to prevent file and line
963 throw TimeSteppingBreakdown{msg};
964 }
965}
966
967template<class TypeTag>
968template<class Solver>
969void
971chopTimeStep_(const double new_time_step)
972{
973 setTimeStep_(new_time_step);
974 if (solverVerbose_()) {
975 const auto msg = fmt::format(fmt::runtime("{}\nTimestep chopped to {} days\n"),
976 this->cause_of_failure_,
977 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
978 OpmLog::problem(msg);
979 }
980}
981
982template<class TypeTag>
983template<class Solver>
984bool
986chopTimeStepOrCloseFailingWells_(const double new_time_step)
987{
988 bool wells_shut = false;
989 // We are below the threshold, and will check if there are any
990 // wells that fails repeatedly (that means that it fails in the last three steps)
991 // we should close rather than chopping again.
992 // If we already have chopped the timestep two times that is
993 // new_time_step < minTimeStepBeforeClosingWells_()*restartFactor_()*restartFactor_()
994 // We also shut wells that fails only on this step.
995 const bool requireRepeatedFailures =
996 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
997 const std::set<std::string> failing_wells =
998 detail::consistentlyFailingWells(solver_().model().stepReports(), requireRepeatedFailures);
999
1000 if (failing_wells.empty()) {
1001 // Found no wells to close, chop the timestep
1002 chopTimeStep_(new_time_step);
1003 } else {
1004 // Close all consistently failing wells that are not under group control
1005 std::vector<std::string> shut_wells;
1006 for (const auto& well : failing_wells) {
1007 const bool was_shut =
1008 solver_().model().wellModel().forceShutWellByName(well,
1009 this->substep_timer_.simulationTimeElapsed(),
1010 /*dont_shut_grup_wells =*/ true);
1011 if (was_shut) {
1012 shut_wells.push_back(well);
1013 }
1014 }
1015 // If no wells are closed we also try to shut wells under group control
1016 if (shut_wells.empty()) {
1017 for (const auto& well : failing_wells) {
1018 const bool was_shut =
1019 solver_().model().wellModel().forceShutWellByName(well,
1020 this->substep_timer_.simulationTimeElapsed(),
1021 /*dont_shut_grup_wells =*/ false);
1022 if (was_shut) {
1023 shut_wells.push_back(well);
1024 }
1025 }
1026 }
1027 // If still no wells are closed we must fall back to chopping again
1028 if (shut_wells.empty()) {
1029 chopTimeStep_(new_time_step);
1030 } else {
1031 wells_shut = true;
1032 if (solverVerbose_()) {
1033 const std::string msg =
1034 fmt::format(fmt::runtime("\nProblematic well(s) were shut: {}"
1035 "(retrying timestep)\n"),
1036 fmt::join(shut_wells, " "));
1037 OpmLog::problem(msg);
1038 }
1039 }
1040 }
1041 return wells_shut;
1042}
1043
1044template<class TypeTag>
1045template<class Solver>
1046boost::posix_time::ptime
1048currentDateTime_() const
1049{
1050 return simulatorTimer_().currentDateTime();
1051}
1052
1053template<class TypeTag>
1054template<class Solver>
1055int
1057getNumIterations_(const SimulatorReportSingle &substep_report) const
1058{
1059 if (useNewtonIteration_()) {
1060 return substep_report.total_newton_iterations;
1061 }
1062 else {
1063 return substep_report.total_linear_iterations;
1064 }
1065}
1066
1067template<class TypeTag>
1068template<class Solver>
1069double
1071growthFactor_() const
1072{
1073 return this->adaptive_time_stepping_.growth_factor_;
1074}
1075
1076template<class TypeTag>
1077template<class Solver>
1078bool
1081{
1082 return adaptive_time_stepping_.ignore_convergence_failure_;
1083}
1084
1085template<class TypeTag>
1086template<class Solver>
1087bool
1090{
1091 return this->substepper_.isReservoirCouplingMaster_();
1092}
1093
1094template<class TypeTag>
1095template<class Solver>
1096bool
1099{
1100 return this->substepper_.isReservoirCouplingSlave_();
1101}
1102
1103template<class TypeTag>
1104template<class Solver>
1105void
1108{
1109#ifdef RESERVOIR_COUPLING_ENABLED
1110 // Clear the first-substep flag after the first runSubStep_() call.
1111 // This ensures that master-slave synchronization only happens once per sync timestep,
1112 // not on retry attempts after convergence-driven timestep chops.
1113 if (isReservoirCouplingMaster_()) {
1114 reservoirCouplingMaster_().setFirstSubstepOfSyncTimestep(false);
1115 }
1116 else if (isReservoirCouplingSlave_()) {
1117 reservoirCouplingSlave_().setFirstSubstepOfSyncTimestep(false);
1118 }
1119#endif
1120 return;
1121}
1122
1123template<class TypeTag>
1124template<class Solver>
1125double
1127maxGrowth_() const
1128{
1129 return this->adaptive_time_stepping_.max_growth_;
1130}
1131
1132template<class TypeTag>
1133template<class Solver>
1134void
1136maybeReportSubStep_(SimulatorReportSingle substep_report) const
1137{
1138 if (timeStepVerbose_()) {
1139 std::ostringstream ss;
1140 substep_report.reportStep(ss);
1141 OpmLog::info(ss.str());
1142 }
1143}
1144
1145template<class TypeTag>
1146template<class Solver>
1147double
1149maybeRestrictTimeStepGrowth_(const double dt, double dt_estimate, const int restarts) const
1150{
1151 // limit the growth of the timestep size by the growth factor
1152 dt_estimate = std::min(dt_estimate, double(maxGrowth_() * dt));
1153 assert(dt_estimate > 0);
1154 // further restrict time step size growth after convergence problems
1155 if (restarts > 0) {
1156 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1157 }
1158
1159 return dt_estimate;
1160}
1161
1162
1163template<class TypeTag>
1164template<class Solver>
1165void
1167maybeUpdateLastSubstepOfSyncTimestep_([[maybe_unused]] const double dt)
1168{
1169#ifdef RESERVOIR_COUPLING_ENABLED
1170 // For reservoir coupling slaves: predict if this substep will complete
1171 // the sync timestep. If so, timeStepSucceeded() will send production
1172 // data to the master (which is blocking on receive after its first substep).
1173 // This is used for summary data synchronization between slaves and master.
1174 if (isReservoirCouplingSlave_()) {
1176 this->substep_timer_.simulationTimeElapsed() + dt,
1177 this->substep_timer_.totalTime()
1178 );
1179 reservoirCouplingSlave_().setLastSubstepOfSyncTimestep(is_last);
1180 }
1181#endif
1182}
1183
1184// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
1185// It has to be called for each substep since TUNING might have been changed for next sub step due
1186// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
1187template<class TypeTag>
1188template<class Solver>
1189void
1192{
1193 // TODO: This function is currently only called if NEXTSTEP is activated from ACTIONX or
1194 // if the WCYCLE keyword needs to modify the current timestep. So this method should rather
1195 // be named maybeUpdateTimeStep_() or similar, since it should not update the tuning. However,
1196 // the current definition of the maybeUpdateTuning_() callback is actually calling
1197 // adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning) which is updating the tuning
1198 // see SimulatorFullyImplicitBlackoil::runStep() for more details.
1199 const auto old_value = suggestedNextTimestep_();
1200 if (this->substepper_.maybeUpdateTuning_(this->substep_timer_.simulationTimeElapsed(),
1201 this->substep_timer_.currentStepLength(),
1202 this->substep_timer_.currentStepNum()))
1203 {
1204 // Either NEXTSTEP and WCYCLE wants to change the current time step, but they cannot
1205 // change the current time step directly. Instead, they change the suggested next time step
1206 // by calling updateNEXTSTEP() via the maybeUpdateTuning() callback. We now need to update
1207 // the current time step to the new suggested time step and reset the suggested time step
1208 // to the old value.
1209 setTimeStep_(suggestedNextTimestep_());
1210 setSuggestedNextStep_(old_value);
1211 }
1212}
1213
1214template<class TypeTag>
1215template<class Solver>
1216double
1219{
1220 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1221}
1222
1223template<class TypeTag>
1224template<class Solver>
1225double
1227minTimeStep_() const
1228{
1229 return this->adaptive_time_stepping_.min_time_step_;
1230}
1231
1232#ifdef RESERVOIR_COUPLING_ENABLED
1233template<class TypeTag>
1234template<class Solver>
1238{
1239 return this->substepper_.reservoirCouplingMaster_();
1240}
1241
1242template<class TypeTag>
1243template<class Solver>
1247{
1248 return this->substepper_.reservoirCouplingSlave_();
1249}
1250#endif
1251
1252template<class TypeTag>
1253template<class Solver>
1254double
1256restartFactor_() const
1257{
1258 return this->adaptive_time_stepping_.restart_factor_;
1259}
1260
1261template<class TypeTag>
1262template<class Solver>
1266{
1267 SimulatorReportSingle substep_report;
1268
1269 auto handleFailure = [this, &substep_report]
1270 (const std::string& failure_reason, const std::exception& e, bool log_exception = true)
1271 {
1272 substep_report = solver_().failureReport();
1273 this->cause_of_failure_ = failure_reason;
1274 if (log_exception && solverVerbose_()) {
1275 OpmLog::debug(std::string("Caught Exception: ") + e.what());
1276 }
1277 };
1278
1279 try {
1280 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1281 if (solverVerbose_()) {
1282 // report number of linear iterations
1283 OpmLog::debug("Overall linear iterations used: "
1284 + std::to_string(substep_report.total_linear_iterations));
1285 }
1286 }
1287 catch (const TooManyIterations& e) {
1288 handleFailure("Solver convergence failure - Iteration limit reached", e);
1289 }
1290 catch (const TimeSteppingBreakdown& e) {
1291 handleFailure(e.what(), e);
1292 }
1293 catch (const ConvergenceMonitorFailure& e) {
1294 handleFailure("Convergence monitor failure", e, /*log_exception=*/false);
1295 }
1296 catch (const LinearSolverProblem& e) {
1297 handleFailure("Linear solver convergence failure", e);
1298 }
1299 catch (const NumericalProblem& e) {
1300 handleFailure("Solver convergence failure - Numerical problem encountered", e);
1301 }
1302 catch (const std::runtime_error& e) {
1303 handleFailure("Runtime error encountered", e);
1304 }
1305 catch (const Dune::ISTLError& e) {
1306 handleFailure("ISTL error - Time step too large", e);
1307 }
1308 catch (const Dune::MatrixBlockError& e) {
1309 handleFailure("Matrix block error", e);
1310 }
1311
1312 return substep_report;
1313}
1314
1315template<class TypeTag>
1316template<class Solver>
1317void
1319setTimeStep_(double dt_estimate)
1320{
1321 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1322}
1323
1324template<class TypeTag>
1325template<class Solver>
1326Solver&
1328solver_() const
1329{
1330 return this->substepper_.solver_;
1331}
1332
1333
1334template<class TypeTag>
1335template<class Solver>
1336int
1338solverRestartMax_() const
1339{
1340 return this->adaptive_time_stepping_.solver_restart_max_;
1341}
1342
1343template<class TypeTag>
1344template<class Solver>
1345void
1347setSuggestedNextStep_(double step)
1348{
1349 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1350}
1351
1352template <class TypeTag>
1353template <class Solver>
1354const SimulatorTimer&
1356simulatorTimer_() const
1357{
1358 return this->substepper_.simulator_timer_;
1359}
1360
1361template <class TypeTag>
1362template <class Solver>
1363bool
1365solverVerbose_() const
1366{
1367 return this->adaptive_time_stepping_.solver_verbose_;
1368}
1369
1370template<class TypeTag>
1371template<class Solver>
1372boost::posix_time::ptime
1374startDateTime_() const
1375{
1376 return simulatorTimer_().startDateTime();
1377}
1378
1379template <class TypeTag>
1380template <class Solver>
1381double
1384{
1385 return this->adaptive_time_stepping_.suggestedNextStep();
1386}
1387
1388template <class TypeTag>
1389template <class Solver>
1390double
1392timeStepControlComputeEstimate_(const double dt, const int iterations,
1393 const AdaptiveSimulatorTimer& substepTimer) const
1394{
1395 // create object to compute the time error, simply forwards the call to the model
1396 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1397 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1398 dt, iterations, relative_change, substepTimer);
1399}
1400
1401template <class TypeTag>
1402template <class Solver>
1403bool
1405timeStepVerbose_() const
1406{
1407 return this->adaptive_time_stepping_.timestep_verbose_;
1408}
1409
1410// The suggested time step is the stepsize that will be used as a first try for
1411// the next sub step. It is updated at the end of each substep. It can also be
1412// updated by the TUNING or NEXTSTEP keywords at the beginning of each report step or
1413// at the beginning of each substep by the ACTIONX keyword (via NEXTSTEP), this is
1414// done by the maybeUpdateTuning_() method which is called at the beginning of each substep
1415// (and the begginning of each report step). Note that the WCYCLE keyword can also update the
1416// suggested time step via the maybeUpdateTuning_() method.
1417template <class TypeTag>
1418template <class Solver>
1419void
1422{
1423 auto suggested_next_step = this->substep_timer_.currentStepLength();
1424 if (! std::isfinite(suggested_next_step)) { // check for NaN
1425 suggested_next_step = this->original_time_step_;
1426 }
1427 if (timeStepVerbose_()) {
1428 std::ostringstream ss;
1429 this->substep_timer_.report(ss);
1430 ss << "Suggested next step size = "
1431 << unit::convert::to(suggested_next_step, unit::day) << " (days)" << std::endl;
1432 OpmLog::debug(ss.str());
1433 }
1434 setSuggestedNextStep_(suggested_next_step);
1435}
1436
1437template <class TypeTag>
1438template <class Solver>
1439bool
1441useNewtonIteration_() const
1442{
1443 return this->adaptive_time_stepping_.use_newton_iteration_;
1444}
1445
1446template <class TypeTag>
1447template <class Solver>
1448double
1450writeOutput_() const
1451{
1452 time::StopWatch perf_timer;
1453 perf_timer.start();
1454 auto& problem = solver_().model().simulator().problem();
1455 problem.writeOutput(true);
1456 return perf_timer.secsSinceStart();
1457}
1458
1459/************************************************
1460 * Private class SolutionTimeErrorSolverWrapper
1461 * **********************************************/
1462
1463template<class TypeTag>
1464template<class Solver>
1467SolutionTimeErrorSolverWrapper(const Solver& solver)
1468 : solver_{solver}
1469{}
1470
1471template<class TypeTag>
1472template<class Solver>
1473double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange() const
1474{
1475 // returns: || u^n+1 - u^n || / || u^n+1 ||
1476 return solver_.model().relativeChange();
1477}
1478
1479} // namespace Opm
1480
1481#endif // OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
Simulation timer for adaptive time stepping.
Definition AdaptiveSimulatorTimer.hpp:41
Definition AdaptiveTimeStepping.hpp:76
double max_growth_
factor that limits the maximum growth of a time step
Definition AdaptiveTimeStepping.hpp:258
double max_time_step_
maximal allowed time step size in days
Definition AdaptiveTimeStepping.hpp:259
bool solver_verbose_
solver verbosity
Definition AdaptiveTimeStepping.hpp:263
int solver_restart_max_
how many restart of solver are allowed
Definition AdaptiveTimeStepping.hpp:262
double timestep_after_event_
suggested size of timestep after an event
Definition AdaptiveTimeStepping.hpp:267
bool ignore_convergence_failure_
continue instead of stop when minimum time step is reached
Definition AdaptiveTimeStepping.hpp:261
double suggested_next_timestep_
suggested size of next timestep
Definition AdaptiveTimeStepping.hpp:265
TimeStepControlType time_step_control_type_
type of time step control object
Definition AdaptiveTimeStepping.hpp:254
bool full_timestep_initially_
beginning with the size of the time step from data file
Definition AdaptiveTimeStepping.hpp:266
SimulatorReport step(const SimulatorTimer &simulator_timer, Solver &solver, const bool is_event, const TuningUpdateCallback &tuning_updater)
step method that acts like the solver::step method in a sub cycle of time steps
Definition AdaptiveTimeStepping_impl.hpp:203
double growth_factor_
factor to multiply time step when solver recovered from failed convergence
Definition AdaptiveTimeStepping.hpp:257
double restart_factor_
factor to multiply time step with when solver fails to converge
Definition AdaptiveTimeStepping.hpp:256
double min_time_step_
minimal allowed time step size before throwing
Definition AdaptiveTimeStepping.hpp:260
TimeStepController time_step_control_
time step control object
Definition AdaptiveTimeStepping.hpp:255
double min_time_step_before_shutting_problematic_wells_
< shut problematic wells when time step size in days are less than this
Definition AdaptiveTimeStepping.hpp:271
bool timestep_verbose_
timestep verbosity
Definition AdaptiveTimeStepping.hpp:264
bool use_newton_iteration_
use newton iteration count for adaptive time step control
Definition AdaptiveTimeStepping.hpp:268
Definition ReservoirCouplingMaster.hpp:38
Definition ReservoirCouplingSlave.hpp:40
Definition SimulatorTimer.hpp:39
TimeStepControlInterface.
Definition TimeStepControlInterface.hpp:51
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
This file provides the infrastructure to retrieve run-time parameters.
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187
static bool compare_gt_or_eq(double a, double b)
Determines if a is greater than b within the specified tolerance.
Definition ReservoirCoupling.cpp:179
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34
Definition SimulatorReport.hpp:122