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/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
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>
66AdaptiveTimeStepping(const UnitSystem& unit_system,
67 const SimulatorReport& report,
68 const double max_next_tstep,
69 const bool terminal_output
70)
71 : time_step_control_{}
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
75 , max_time_step_{
76 Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60} // 365.25
77 , min_time_step_{
78 unit_system.to_si(UnitSystem::measure::time,
79 Parameters::Get<Parameters::SolverMinTimeStep<Scalar>>())} // 1e-12;
80 , ignore_convergence_failure_{
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
85 , suggested_next_timestep_{
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
89 , timestep_after_event_{
90 Parameters::Get<Parameters::TimeStepAfterEventInDays<Scalar>>() * 24 * 60 * 60} // 1e30
91 , use_newton_iteration_{false}
92 , min_time_step_before_shutting_problematic_wells_{
93 Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day}
94 , report_(report)
95{
96 init_(unit_system);
97}
98
105template<class 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)
113 : time_step_control_{}
114 , restart_factor_{tuning.TSFCNV}
115 , growth_factor_{tuning.TFDIFF}
116 , max_growth_{tuning.TSFMAX}
117 , max_time_step_{tuning.TSMAXZ} // 365.25
118 , min_time_step_{tuning.TSMINZ} // 0.1;
119 , ignore_convergence_failure_{true}
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
123 , suggested_next_timestep_{
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}
129 , min_time_step_before_shutting_problematic_wells_{
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_) {
150 result = castAndComp<HardcodedTimeStepControl>(rhs);
151 break;
153 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
154 break;
156 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
157 break;
159 result = castAndComp<PIDTimeStepControl>(rhs);
160 break;
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>();
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_) {
223 allocAndSerialize<HardcodedTimeStepControl>(serializer);
224 break;
226 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
227 break;
229 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
230 break;
232 allocAndSerialize<PIDTimeStepControl>(serializer);
233 break;
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
372AdaptiveTimeStepping<TypeTag>::
373castAndComp(const AdaptiveTimeStepping<TypeTag>& Rhs) const
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
382AdaptiveTimeStepping<TypeTag>::
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>
403AdaptiveTimeStepping<TypeTag>
404AdaptiveTimeStepping<TypeTag>::
405serializationTestObject_()
406{
407 AdaptiveTimeStepping<TypeTag> result;
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>
471AdaptiveTimeStepping<TypeTag>&
472AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
473getAdaptiveTimerStepper()
474{
475 return adaptive_time_stepping_;
476}
477
478template<class TypeTag>
479template<class Solver>
480SimulatorReport
481AdaptiveTimeStepping<TypeTag>::SubStepper<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
507AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
508isReservoirCouplingMaster_() const
509{
510 return this->solver_.model().simulator().reservoirCouplingMaster() != nullptr;
511}
512
513template<class TypeTag>
514template<class Solver>
515bool
516AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
517isReservoirCouplingSlave_() const
518{
519 return this->solver_.model().simulator().reservoirCouplingSlave() != nullptr;
520}
521#endif
522
523template<class TypeTag>
524template<class Solver>
525void
526AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
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
540AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
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
549AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
550maxTimeStep_() const
551{
552 return this->adaptive_time_stepping_.max_time_step_;
553}
554
555template <class TypeTag>
556template <class Solver>
557SimulatorReport
558AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
559runStepOriginal_()
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>
582ReservoirCouplingMaster<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
583AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
584reservoirCouplingMaster_()
585{
586 return *(this->solver_.model().simulator().reservoirCouplingMaster());
587}
588#endif
589
590#ifdef RESERVOIR_COUPLING_ENABLED
591template <class TypeTag>
592template <class Solver>
593ReservoirCouplingSlave<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
594AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
595reservoirCouplingSlave_()
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>
634SimulatorReport
635AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
636runStepReservoirCouplingMaster_()
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_().receiveNextReportDateFromSlaves();
652 bool start_of_report_step = (iteration == 0);
653 if (start_of_report_step) {
654 reservoirCouplingMaster_().initStartOfReportStep(report_step_idx);
655 }
656 current_step_length = reservoirCouplingMaster_().maybeChopSubStep(
657 current_step_length, current_time);
658 reservoirCouplingMaster_().sendNextTimeStepToSlaves(current_step_length);
659 if (start_of_report_step) {
660 maybeUpdateTuning_(current_time, current_step_length, /*substep=*/0);
661 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(current_step_length);
662 }
663 AdaptiveSimulatorTimer substep_timer{
664 this->simulator_timer_.startDateTime(),
665 /*stepLength=*/current_step_length,
666 /*elapsedTime=*/current_time,
667 /*timeStepEstimate=*/suggestedNextTimestep_(),
668 /*reportStep=*/this->simulator_timer_.reportStepNum(),
669 maxTimeStep_()
670 };
671 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
672 current_time + current_step_length, step_end_time
673 );
674 // Mark this as the first substep of the "sync" timestep. This flag controls
675 // whether master-slave data exchange should occur in beginTimeStep() in the well model.
676 // It will be cleared after the first runSubStep_() call.
677 reservoirCouplingMaster_().setFirstSubstepOfSyncTimestep(true);
678 SubStepIteration<Solver> substepIteration{*this, substep_timer, current_step_length, final_step};
679 const auto sub_steps_report = substepIteration.run();
680 report += sub_steps_report;
681 current_time += current_step_length;
682 if (final_step) {
683 break;
684 }
685 iteration++;
686 }
687 return report;
688}
689#endif
690
691#ifdef RESERVOIR_COUPLING_ENABLED
692template <class TypeTag>
693template <class Solver>
694SimulatorReport
695AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
696runStepReservoirCouplingSlave_()
697{
698 int iteration = 0;
699 const double original_time_step = this->simulator_timer_.currentStepLength();
700 double current_time{this->simulator_timer_.simulationTimeElapsed()};
701 double step_end_time = current_time + original_time_step;
702 SimulatorReport report;
703 auto report_step_idx = this->simulator_timer_.currentStepNum();
704 if (report_step_idx == 0 && iteration == 0) {
705 reservoirCouplingSlave_().initTimeStepping();
706 }
707 while (true) {
708 bool start_of_report_step = (iteration == 0);
709 reservoirCouplingSlave_().sendNextReportDateToMasterProcess();
710 const auto timestep = reservoirCouplingSlave_().receiveNextTimeStepFromMaster();
711 if (start_of_report_step) {
712 maybeUpdateTuning_(current_time, original_time_step, /*substep=*/0);
713 maybeModifySuggestedTimeStepAtBeginningOfReportStep_(timestep);
714 }
715 AdaptiveSimulatorTimer substep_timer{
716 this->simulator_timer_.startDateTime(),
717 /*step_length=*/timestep,
718 /*elapsed_time=*/current_time,
719 /*time_step_estimate=*/suggestedNextTimestep_(),
720 this->simulator_timer_.reportStepNum(),
721 maxTimeStep_()
722 };
723 const bool final_step = ReservoirCoupling::Seconds::compare_gt_or_eq(
724 current_time + timestep, step_end_time
725 );
726 // Mark this as the first substep of the "sync" timestep. This flag controls
727 // whether master-slave data exchange should occur in beginTimeStep() in the well model.
728 // It will be cleared after the first runSubStep_() call.
729 reservoirCouplingSlave_().setFirstSubstepOfSyncTimestep(true);
730 SubStepIteration<Solver> substepIteration{*this, substep_timer, timestep, final_step};
731 const auto sub_steps_report = substepIteration.run();
732 report += sub_steps_report;
733 current_time += timestep;
734 if (final_step) {
735 break;
736 }
737 iteration++;
738 }
739 return report;
740}
741
742#endif
743
744template <class TypeTag>
745template <class Solver>
746double
747AdaptiveTimeStepping<TypeTag>::SubStepper<Solver>::
748suggestedNextTimestep_() const
749{
750 return this->adaptive_time_stepping_.suggestedNextStep();
751}
752
753
754
755/************************************************
756 * Private class SubStepIteration public methods
757 ************************************************/
758
759template<class TypeTag>
760template<class Solver>
761AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
762SubStepIteration(SubStepper<Solver>& substepper,
763 AdaptiveSimulatorTimer& substep_timer,
764 const double original_time_step,
765 bool final_step)
766 : substepper_{substepper}
767 , substep_timer_{substep_timer}
768 , original_time_step_{original_time_step}
769 , final_step_{final_step}
770 , adaptive_time_stepping_{substepper.getAdaptiveTimerStepper()}
771{
772}
773
774template <class TypeTag>
775template <class Solver>
776SimulatorReport
777AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
778run()
779{
780 auto& simulator = solver_().model().simulator();
781 auto& problem = simulator.problem();
782 // counter for solver restarts
783 int restarts = 0;
784 SimulatorReport report;
785
786 // sub step time loop
787 while (!this->substep_timer_.done()) {
788 // if we just chopped the timestep due to convergence i.e. restarts>0
789 // we dont want to update the next timestep based on Tuning
790 if (restarts == 0) {
791 maybeUpdateTuningAndTimeStep_();
792 }
793 const double dt = this->substep_timer_.currentStepLength();
794 if (timeStepVerbose_()) {
795 detail::logTimer(this->substep_timer_);
796 }
797
798 auto substep_report = runSubStep_();
799 markFirstSubStepAsFinished_(); // Needed for reservoir coupling
800
801 if (substep_report.converged || checkContinueOnUnconvergedSolution_(dt)) {
802 Dune::Timer perfTimer;
803 perfTimer.start();
804 // Pass substep to eclwriter for summary output
805 problem.setSubStepReport(substep_report);
806 auto& full_report = adaptive_time_stepping_.report();
807 full_report += substep_report;
808 problem.setSimulationReport(full_report);
809 problem.endTimeStep();
810 substep_report.pre_post_time += perfTimer.stop();
811
812 report += substep_report;
813
814 ++this->substep_timer_; // advance by current dt
815
816 const int iterations = getNumIterations_(substep_report);
817 auto dt_estimate = timeStepControlComputeEstimate_(
818 dt, iterations, this->substep_timer_);
819
820 assert(dt_estimate > 0);
821 dt_estimate = maybeRestrictTimeStepGrowth_(dt, dt_estimate, restarts);
822 restarts = 0; // solver converged, reset restarts counter
823
824 maybeReportSubStep_(substep_report);
825 if (this->final_step_ && this->substep_timer_.done()) {
826 // if the time step is done we do not need to write it as this will be done
827 // by the simulator anyway.
828 }
829 else {
830 report.success.output_write_time += writeOutput_();
831 }
832
833 // set new time step length
834 setTimeStep_(dt_estimate);
835
836 report.success.converged = this->substep_timer_.done();
837 this->substep_timer_.setLastStepFailed(false);
838 }
839 else { // in case of no convergence or time step tolerance test failure
840 report += substep_report;
841 this->substep_timer_.setLastStepFailed(true);
842 checkTimeStepMaxRestartLimit_(restarts);
843
844 double new_time_step = restartFactor_() * dt;
845 if (substep_report.time_step_rejected) {
846 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>();
847 const double safetyFactor = Parameters::Get<Parameters::TimeStepControlSafetyFactor>();
848 const double temp_time_step = std::sqrt(safetyFactor * tol / solver_().model().relativeChange()) * dt;
849 if (temp_time_step < dt) { // added in case suggested time step is not a reduction
850 new_time_step = temp_time_step;
851 }
852 }
853 checkTimeStepMinLimit_(new_time_step);
854 bool wells_shut = false;
855 if (new_time_step > minTimeStepBeforeClosingWells_()) {
856 chopTimeStep_(new_time_step);
857 } else {
858 wells_shut = chopTimeStepOrCloseFailingWells_(new_time_step);
859 }
860 if (wells_shut) {
861 setTimeStep_(dt); // retry the old timestep
862 }
863 else {
864 restarts++; // only increase if no wells were shut
865 }
866 }
867 problem.setNextTimeStepSize(this->substep_timer_.currentStepLength());
868 }
869 updateSuggestedNextStep_();
870 return report;
871}
872
873
874/************************************************
875 * Private class SubStepIteration private methods
876 ************************************************/
877
878
879template<class TypeTag>
880template<class Solver>
881bool
882AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
883checkContinueOnUnconvergedSolution_(double dt) const
884{
885 const bool continue_on_uncoverged_solution = ignoreConvergenceFailure_() && dt <= minTimeStep_();
886 if (continue_on_uncoverged_solution && solverVerbose_()) {
887 // NOTE: This method is only called if the solver failed to converge.
888 const auto msg = fmt::format(
889 "Solver failed to converge but timestep {} is smaller or equal to {}\n"
890 "which is the minimum threshold given by option --solver-min-time-step\n",
891 dt, minTimeStep_()
892 );
893 OpmLog::problem(msg);
894 }
895 return continue_on_uncoverged_solution;
896}
897
898template<class TypeTag>
899template<class Solver>
900void
901AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
902checkTimeStepMaxRestartLimit_(const int restarts) const
903{
904 // If we have restarted (i.e. cut the timestep) too
905 // many times, we have failed and throw an exception.
906 if (restarts >= solverRestartMax_()) {
907 const auto msg = fmt::format(
908 "Solver failed to converge after cutting timestep {} times.", restarts
909 );
910 if (solverVerbose_()) {
911 OpmLog::error(msg);
912 }
913 // Use throw directly to prevent file and line
914 throw TimeSteppingBreakdown{msg};
915 }
916}
917
918template<class TypeTag>
919template<class Solver>
920void
921AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
922checkTimeStepMinLimit_(const double new_time_step) const
923{
924 using Meas = UnitSystem::measure;
925 // If we have restarted (i.e. cut the timestep) too
926 // much, we have failed and throw an exception.
927 if (new_time_step < minTimeStep_()) {
928 auto msg = fmt::format("Solver failed to converge after cutting timestep to ");
929 if (Parameters::Get<Parameters::EnableTuning>()) {
930 const UnitSystem& unit_system = solver_().model().simulator().vanguard().eclState().getDeckUnitSystem();
931 msg += fmt::format(
932 "{:.3E} {}\nwhich is the minimum threshold given by the TUNING keyword\n",
933 unit_system.from_si(Meas::time, minTimeStep_()),
934 unit_system.name(Meas::time)
935 );
936 }
937 else {
938 msg += fmt::format(
939 "{:.3E} DAYS\nwhich is the minimum threshold given by option --solver-min-time-step\n",
940 minTimeStep_() / 86400.0
941 );
942 }
943 if (solverVerbose_()) {
944 OpmLog::error(msg);
945 }
946 // Use throw directly to prevent file and line
947 throw TimeSteppingBreakdown{msg};
948 }
949}
950
951template<class TypeTag>
952template<class Solver>
953void
954AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
955chopTimeStep_(const double new_time_step)
956{
957 setTimeStep_(new_time_step);
958 if (solverVerbose_()) {
959 const auto msg = fmt::format("{}\nTimestep chopped to {} days\n",
960 this->cause_of_failure_,
961 unit::convert::to(this->substep_timer_.currentStepLength(), unit::day));
962 OpmLog::problem(msg);
963 }
964}
965
966template<class TypeTag>
967template<class Solver>
968bool
969AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
970chopTimeStepOrCloseFailingWells_(const double new_time_step)
971{
972 bool wells_shut = false;
973 // We are below the threshold, and will check if there are any
974 // wells that fails repeatedly (that means that it fails in the last three steps)
975 // we should close rather than chopping again.
976 // If we already have chopped the timestep two times that is
977 // new_time_step < minTimeStepBeforeClosingWells_()*restartFactor_()*restartFactor_()
978 // We also shut wells that fails only on this step.
979 const bool requireRepeatedFailures =
980 new_time_step > (minTimeStepBeforeClosingWells_() * restartFactor_() * restartFactor_());
981 const std::set<std::string> failing_wells =
982 detail::consistentlyFailingWells(solver_().model().stepReports(), requireRepeatedFailures);
983
984 if (failing_wells.empty()) {
985 // Found no wells to close, chop the timestep
986 chopTimeStep_(new_time_step);
987 } else {
988 // Close all consistently failing wells that are not under group control
989 std::vector<std::string> shut_wells;
990 for (const auto& well : failing_wells) {
991 const bool was_shut =
992 solver_().model().wellModel().forceShutWellByName(well,
993 this->substep_timer_.simulationTimeElapsed(),
994 /*dont_shut_grup_wells =*/ true);
995 if (was_shut) {
996 shut_wells.push_back(well);
997 }
998 }
999 // If no wells are closed we also try to shut wells under group control
1000 if (shut_wells.empty()) {
1001 for (const auto& well : failing_wells) {
1002 const bool was_shut =
1003 solver_().model().wellModel().forceShutWellByName(well,
1004 this->substep_timer_.simulationTimeElapsed(),
1005 /*dont_shut_grup_wells =*/ false);
1006 if (was_shut) {
1007 shut_wells.push_back(well);
1008 }
1009 }
1010 }
1011 // If still no wells are closed we must fall back to chopping again
1012 if (shut_wells.empty()) {
1013 chopTimeStep_(new_time_step);
1014 } else {
1015 wells_shut = true;
1016 if (solverVerbose_()) {
1017 const std::string msg =
1018 fmt::format("\nProblematic well(s) were shut: {}"
1019 "(retrying timestep)\n",
1020 fmt::join(shut_wells, " "));
1021 OpmLog::problem(msg);
1022 }
1023 }
1024 }
1025 return wells_shut;
1026}
1027
1028template<class TypeTag>
1029template<class Solver>
1030boost::posix_time::ptime
1031AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1032currentDateTime_() const
1033{
1034 return simulatorTimer_().currentDateTime();
1035}
1036
1037template<class TypeTag>
1038template<class Solver>
1039int
1040AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1041getNumIterations_(const SimulatorReportSingle &substep_report) const
1042{
1043 if (useNewtonIteration_()) {
1044 return substep_report.total_newton_iterations;
1045 }
1046 else {
1047 return substep_report.total_linear_iterations;
1048 }
1049}
1050
1051template<class TypeTag>
1052template<class Solver>
1053double
1054AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1055growthFactor_() const
1056{
1057 return this->adaptive_time_stepping_.growth_factor_;
1058}
1059
1060template<class TypeTag>
1061template<class Solver>
1062bool
1063AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1064ignoreConvergenceFailure_() const
1065{
1066 return adaptive_time_stepping_.ignore_convergence_failure_;
1067}
1068
1069template<class TypeTag>
1070template<class Solver>
1071bool
1072AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1073isReservoirCouplingMaster_() const
1074{
1075 return this->substepper_.isReservoirCouplingMaster_();
1076}
1077
1078template<class TypeTag>
1079template<class Solver>
1080bool
1081AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1082isReservoirCouplingSlave_() const
1083{
1084 return this->substepper_.isReservoirCouplingSlave_();
1085}
1086
1087template<class TypeTag>
1088template<class Solver>
1089void
1090AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1091markFirstSubStepAsFinished_() const
1092{
1093#ifdef RESERVOIR_COUPLING_ENABLED
1094 // Clear the first-substep flag after the first runSubStep_() call.
1095 // This ensures that master-slave synchronization only happens once per sync timestep,
1096 // not on retry attempts after convergence-driven timestep chops.
1097 if (isReservoirCouplingMaster_()) {
1098 reservoirCouplingMaster_().setFirstSubstepOfSyncTimestep(false);
1099 }
1100 else if (isReservoirCouplingSlave_()) {
1101 reservoirCouplingSlave_().setFirstSubstepOfSyncTimestep(false);
1102 }
1103#endif
1104 return;
1105}
1106
1107template<class TypeTag>
1108template<class Solver>
1109double
1110AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1111maxGrowth_() const
1112{
1113 return this->adaptive_time_stepping_.max_growth_;
1114}
1115
1116template<class TypeTag>
1117template<class Solver>
1118void
1119AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1120maybeReportSubStep_(SimulatorReportSingle substep_report) const
1121{
1122 if (timeStepVerbose_()) {
1123 std::ostringstream ss;
1124 substep_report.reportStep(ss);
1125 OpmLog::info(ss.str());
1126 }
1127}
1128
1129template<class TypeTag>
1130template<class Solver>
1131double
1132AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1133maybeRestrictTimeStepGrowth_(const double dt, double dt_estimate, const int restarts) const
1134{
1135 // limit the growth of the timestep size by the growth factor
1136 dt_estimate = std::min(dt_estimate, double(maxGrowth_() * dt));
1137 assert(dt_estimate > 0);
1138 // further restrict time step size growth after convergence problems
1139 if (restarts > 0) {
1140 dt_estimate = std::min(growthFactor_() * dt, dt_estimate);
1141 }
1142
1143 return dt_estimate;
1144}
1145
1146// The maybeUpdateTuning_() lambda callback is defined in SimulatorFullyImplicitBlackoil::runStep()
1147// It has to be called for each substep since TUNING might have been changed for next sub step due
1148// to ACTIONX (via NEXTSTEP) or WCYCLE keywords.
1149template<class TypeTag>
1150template<class Solver>
1151void
1152AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1153maybeUpdateTuningAndTimeStep_()
1154{
1155 // TODO: This function is currently only called if NEXTSTEP is activated from ACTIONX or
1156 // if the WCYCLE keyword needs to modify the current timestep. So this method should rather
1157 // be named maybeUpdateTimeStep_() or similar, since it should not update the tuning. However,
1158 // the current definition of the maybeUpdateTuning_() callback is actually calling
1159 // adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning) which is updating the tuning
1160 // see SimulatorFullyImplicitBlackoil::runStep() for more details.
1161 const auto old_value = suggestedNextTimestep_();
1162 if (this->substepper_.maybeUpdateTuning_(this->substep_timer_.simulationTimeElapsed(),
1163 this->substep_timer_.currentStepLength(),
1164 this->substep_timer_.currentStepNum()))
1165 {
1166 // Either NEXTSTEP and WCYCLE wants to change the current time step, but they cannot
1167 // change the current time step directly. Instead, they change the suggested next time step
1168 // by calling updateNEXTSTEP() via the maybeUpdateTuning() callback. We now need to update
1169 // the current time step to the new suggested time step and reset the suggested time step
1170 // to the old value.
1171 setTimeStep_(suggestedNextTimestep_());
1172 setSuggestedNextStep_(old_value);
1173 }
1174}
1175
1176template<class TypeTag>
1177template<class Solver>
1178double
1179AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1180minTimeStepBeforeClosingWells_() const
1181{
1182 return this->adaptive_time_stepping_.min_time_step_before_shutting_problematic_wells_;
1183}
1184
1185template<class TypeTag>
1186template<class Solver>
1187double
1188AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1189minTimeStep_() const
1190{
1191 return this->adaptive_time_stepping_.min_time_step_;
1192}
1193
1194#ifdef RESERVOIR_COUPLING_ENABLED
1195template<class TypeTag>
1196template<class Solver>
1197ReservoirCouplingMaster<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
1198AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1199reservoirCouplingMaster_() const
1200{
1201 return this->substepper_.reservoirCouplingMaster_();
1202}
1203
1204template<class TypeTag>
1205template<class Solver>
1206ReservoirCouplingSlave<typename AdaptiveTimeStepping<TypeTag>::Scalar>&
1207AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1208reservoirCouplingSlave_() const
1209{
1210 return this->substepper_.reservoirCouplingSlave_();
1211}
1212#endif
1213
1214template<class TypeTag>
1215template<class Solver>
1216double
1217AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1218restartFactor_() const
1219{
1220 return this->adaptive_time_stepping_.restart_factor_;
1221}
1222
1223template<class TypeTag>
1224template<class Solver>
1225SimulatorReportSingle
1226AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1227runSubStep_()
1228{
1229 SimulatorReportSingle substep_report;
1230
1231 auto handleFailure = [this, &substep_report]
1232 (const std::string& failure_reason, const std::exception& e, bool log_exception = true)
1233 {
1234 substep_report = solver_().failureReport();
1235 this->cause_of_failure_ = failure_reason;
1236 if (log_exception && solverVerbose_()) {
1237 OpmLog::debug(std::string("Caught Exception: ") + e.what());
1238 }
1239 };
1240
1241 try {
1242 substep_report = solver_().step(this->substep_timer_, &this->adaptive_time_stepping_.timeStepControl());
1243 if (solverVerbose_()) {
1244 // report number of linear iterations
1245 OpmLog::debug("Overall linear iterations used: "
1246 + std::to_string(substep_report.total_linear_iterations));
1247 }
1248 }
1249 catch (const TooManyIterations& e) {
1250 handleFailure("Solver convergence failure - Iteration limit reached", e);
1251 }
1252 catch (const TimeSteppingBreakdown& e) {
1253 handleFailure(e.what(), e);
1254 }
1255 catch (const ConvergenceMonitorFailure& e) {
1256 handleFailure("Convergence monitor failure", e, /*log_exception=*/false);
1257 }
1258 catch (const LinearSolverProblem& e) {
1259 handleFailure("Linear solver convergence failure", e);
1260 }
1261 catch (const NumericalProblem& e) {
1262 handleFailure("Solver convergence failure - Numerical problem encountered", e);
1263 }
1264 catch (const std::runtime_error& e) {
1265 handleFailure("Runtime error encountered", e);
1266 }
1267 catch (const Dune::ISTLError& e) {
1268 handleFailure("ISTL error - Time step too large", e);
1269 }
1270 catch (const Dune::MatrixBlockError& e) {
1271 handleFailure("Matrix block error", e);
1272 }
1273
1274 return substep_report;
1275}
1276
1277template<class TypeTag>
1278template<class Solver>
1279void
1280AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1281setTimeStep_(double dt_estimate)
1282{
1283 this->substep_timer_.provideTimeStepEstimate(dt_estimate);
1284}
1285
1286template<class TypeTag>
1287template<class Solver>
1288Solver&
1289AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1290solver_() const
1291{
1292 return this->substepper_.solver_;
1293}
1294
1295
1296template<class TypeTag>
1297template<class Solver>
1298int
1299AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1300solverRestartMax_() const
1301{
1302 return this->adaptive_time_stepping_.solver_restart_max_;
1303}
1304
1305template<class TypeTag>
1306template<class Solver>
1307void
1308AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1309setSuggestedNextStep_(double step)
1310{
1311 this->adaptive_time_stepping_.setSuggestedNextStep(step);
1312}
1313
1314template <class TypeTag>
1315template <class Solver>
1316const SimulatorTimer&
1317AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1318simulatorTimer_() const
1319{
1320 return this->substepper_.simulator_timer_;
1321}
1322
1323template <class TypeTag>
1324template <class Solver>
1325bool
1326AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1327solverVerbose_() const
1328{
1329 return this->adaptive_time_stepping_.solver_verbose_;
1330}
1331
1332template<class TypeTag>
1333template<class Solver>
1334boost::posix_time::ptime
1335AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1336startDateTime_() const
1337{
1338 return simulatorTimer_().startDateTime();
1339}
1340
1341template <class TypeTag>
1342template <class Solver>
1343double
1344AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1345suggestedNextTimestep_() const
1346{
1347 return this->adaptive_time_stepping_.suggestedNextStep();
1348}
1349
1350template <class TypeTag>
1351template <class Solver>
1352double
1353AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1354timeStepControlComputeEstimate_(const double dt, const int iterations,
1355 const AdaptiveSimulatorTimer& substepTimer) const
1356{
1357 // create object to compute the time error, simply forwards the call to the model
1358 const SolutionTimeErrorSolverWrapper<Solver> relative_change{solver_()};
1359 return this->adaptive_time_stepping_.time_step_control_->computeTimeStepSize(
1360 dt, iterations, relative_change, substepTimer);
1361}
1362
1363template <class TypeTag>
1364template <class Solver>
1365bool
1366AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1367timeStepVerbose_() const
1368{
1369 return this->adaptive_time_stepping_.timestep_verbose_;
1370}
1371
1372// The suggested time step is the stepsize that will be used as a first try for
1373// the next sub step. It is updated at the end of each substep. It can also be
1374// updated by the TUNING or NEXTSTEP keywords at the beginning of each report step or
1375// at the beginning of each substep by the ACTIONX keyword (via NEXTSTEP), this is
1376// done by the maybeUpdateTuning_() method which is called at the beginning of each substep
1377// (and the begginning of each report step). Note that the WCYCLE keyword can also update the
1378// suggested time step via the maybeUpdateTuning_() method.
1379template <class TypeTag>
1380template <class Solver>
1381void
1382AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1383updateSuggestedNextStep_()
1384{
1385 auto suggested_next_step = this->substep_timer_.currentStepLength();
1386 if (! std::isfinite(suggested_next_step)) { // check for NaN
1387 suggested_next_step = this->original_time_step_;
1388 }
1389 if (timeStepVerbose_()) {
1390 std::ostringstream ss;
1391 this->substep_timer_.report(ss);
1392 ss << "Suggested next step size = "
1393 << unit::convert::to(suggested_next_step, unit::day) << " (days)" << std::endl;
1394 OpmLog::debug(ss.str());
1395 }
1396 setSuggestedNextStep_(suggested_next_step);
1397}
1398
1399template <class TypeTag>
1400template <class Solver>
1401bool
1402AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1403useNewtonIteration_() const
1404{
1405 return this->adaptive_time_stepping_.use_newton_iteration_;
1406}
1407
1408template <class TypeTag>
1409template <class Solver>
1410double
1411AdaptiveTimeStepping<TypeTag>::SubStepIteration<Solver>::
1412writeOutput_() const
1413{
1414 time::StopWatch perf_timer;
1415 perf_timer.start();
1416 auto& problem = solver_().model().simulator().problem();
1417 problem.writeOutput(true);
1418 return perf_timer.secsSinceStart();
1419}
1420
1421/************************************************
1422 * Private class SolutionTimeErrorSolverWrapper
1423 * **********************************************/
1424
1425template<class TypeTag>
1426template<class Solver>
1427AdaptiveTimeStepping<TypeTag>::
1428SolutionTimeErrorSolverWrapper<Solver>::
1429SolutionTimeErrorSolverWrapper(const Solver& solver)
1430 : solver_{solver}
1431{}
1432
1433template<class TypeTag>
1434template<class Solver>
1435double AdaptiveTimeStepping<TypeTag>::SolutionTimeErrorSolverWrapper<Solver>::relativeChange() const
1436{
1437 // returns: || u^n+1 - u^n || / || u^n+1 ||
1438 return solver_.model().relativeChange();
1439}
1440
1441} // namespace Opm
1442
1443#endif // OPM_ADAPTIVE_TIME_STEPPING_IMPL_HPP
Definition: AdaptiveTimeStepping.hpp:76
double max_growth_
factor that limits the maximum growth of a time step
Definition: AdaptiveTimeStepping.hpp:257
double max_time_step_
maximal allowed time step size in days
Definition: AdaptiveTimeStepping.hpp:258
bool solver_verbose_
solver verbosity
Definition: AdaptiveTimeStepping.hpp:262
int solver_restart_max_
how many restart of solver are allowed
Definition: AdaptiveTimeStepping.hpp:261
double timestep_after_event_
suggested size of timestep after an event
Definition: AdaptiveTimeStepping.hpp:266
void init_(const UnitSystem &unitSystem)
Definition: AdaptiveTimeStepping_impl.hpp:435
void setSuggestedNextStep(const double x)
Definition: AdaptiveTimeStepping_impl.hpp:306
double suggestedNextStep() const
Definition: AdaptiveTimeStepping_impl.hpp:314
std::function< bool(const double, const double, const int)> TuningUpdateCallback
Definition: AdaptiveTimeStepping.hpp:78
bool operator==(const AdaptiveTimeStepping< TypeTag > &rhs) const
Definition: AdaptiveTimeStepping_impl.hpp:139
static AdaptiveTimeStepping< TypeTag > serializationTestObjectSimple()
Definition: AdaptiveTimeStepping_impl.hpp:289
bool ignore_convergence_failure_
continue instead of stop when minimum time step is reached
Definition: AdaptiveTimeStepping.hpp:260
void serializeOp(Serializer &serializer)
Definition: AdaptiveTimeStepping_impl.hpp:218
void updateTUNING(double max_next_tstep, const Tuning &tuning)
Definition: AdaptiveTimeStepping_impl.hpp:343
TimeStepControlType time_step_control_type_
type of time step control object
Definition: AdaptiveTimeStepping.hpp:253
const TimeStepControlInterface & timeStepControl() const
Definition: AdaptiveTimeStepping_impl.hpp:322
bool full_timestep_initially_
beginning with the size of the time step from data file
Definition: AdaptiveTimeStepping.hpp:265
SimulatorReport step(const SimulatorTimer &simulator_timer, Solver &solver, const bool is_event, const TuningUpdateCallback &tuning_updater)
step method that acts like the solver::step method in a sub cycle of time steps
Definition: AdaptiveTimeStepping_impl.hpp:203
double growth_factor_
factor to multiply time step when solver recovered from failed convergence
Definition: AdaptiveTimeStepping.hpp:256
double restart_factor_
factor to multiply time step with when solver fails to converge
Definition: AdaptiveTimeStepping.hpp:255
double min_time_step_
minimal allowed time step size before throwing
Definition: AdaptiveTimeStepping.hpp:259
void updateNEXTSTEP(double max_next_tstep)
Definition: AdaptiveTimeStepping_impl.hpp:331
static AdaptiveTimeStepping< TypeTag > serializationTestObjectHardcoded()
Definition: AdaptiveTimeStepping_impl.hpp:265
TimeStepController time_step_control_
time step control object
Definition: AdaptiveTimeStepping.hpp:254
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPIDIt()
Definition: AdaptiveTimeStepping_impl.hpp:281
double min_time_step_before_shutting_problematic_wells_
< shut problematic wells when time step size in days are less than this
Definition: AdaptiveTimeStepping.hpp:270
static AdaptiveTimeStepping< TypeTag > serializationTestObject3rdOrder()
Definition: AdaptiveTimeStepping_impl.hpp:297
SimulatorReport & report()
Definition: AdaptiveTimeStepping_impl.hpp:257
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPID()
Definition: AdaptiveTimeStepping_impl.hpp:273
static void registerParameters()
Definition: AdaptiveTimeStepping_impl.hpp:185
bool use_newton_iteration_
use newton iteration count for adaptive time step control
Definition: AdaptiveTimeStepping.hpp:267
Definition: SimulatorTimer.hpp:39
Definition: TimeStepControlInterface.hpp:51
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:187
void logTimer(const AdaptiveSimulatorTimer &substep_timer)
void registerAdaptiveParameters()
std::set< std::string > consistentlyFailingWells(const std::vector< StepReport > &sr, bool requireRepeatedFailures)
std::tuple< TimeStepControlType, std::unique_ptr< TimeStepControlInterface >, bool > createController(const UnitSystem &unitSystem)
Definition: blackoilbioeffectsmodules.hh:43
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
This file provides the infrastructure to retrieve run-time parameters.
static bool compare_gt_or_eq(double a, double b)
Determines if a is greater than b within the specified tolerance.
Definition: SimulatorReport.hpp:122