AdaptiveTimeStepping.hpp
Go to the documentation of this file.
1/*
2*/
3#ifndef OPM_ADAPTIVE_TIME_STEPPING_HPP
4#define OPM_ADAPTIVE_TIME_STEPPING_HPP
5
6#include <dune/common/version.hh>
7#if DUNE_VERSION_NEWER(DUNE_ISTL, 2, 8)
8#include <dune/istl/istlexception.hh>
9#else
10#include <dune/istl/ilu.hh>
11#endif
12
13#include <opm/common/Exceptions.hpp>
14#include <opm/common/ErrorMacros.hpp>
15#include <opm/common/OpmLog/OpmLog.hpp>
16
17#include <opm/grid/utility/StopWatch.hpp>
18
19#include <opm/input/eclipse/Units/Units.hpp>
20#include <opm/input/eclipse/Units/UnitSystem.hpp>
21
22#include <opm/input/eclipse/Schedule/Tuning.hpp>
23
27
34
36
37#include <fmt/format.h>
38
39#include <algorithm>
40#include <cassert>
41#include <cmath>
42#include <functional>
43#include <memory>
44#include <set>
45#include <sstream>
46#include <stdexcept>
47#include <string>
48#include <vector>
49
50namespace Opm::Parameters {
51
52struct SolverContinueOnConvergenceFailure { static constexpr bool value = false; };
53struct SolverMaxRestarts { static constexpr int value = 10; };
54struct SolverVerbosity { static constexpr int value = 1; };
55struct TimeStepVerbosity { static constexpr int value = 1; };
56struct InitialTimeStepInDays { static constexpr double value = 1.0; };
57struct FullTimeStepInitially { static constexpr bool value = false; };
58struct TimeStepControl { static constexpr auto value = "pid+newtoniteration"; };
59struct TimeStepControlTolerance { static constexpr double value = 1e-1; };
60struct TimeStepControlTargetIterations { static constexpr int value = 30; };
61struct TimeStepControlTargetNewtonIterations { static constexpr int value = 8; };
62struct TimeStepControlDecayRate { static constexpr double value = 0.75; };
63struct TimeStepControlGrowthRate { static constexpr double value = 1.25; };
64struct TimeStepControlDecayDampingFactor { static constexpr double value = 1.0; };
65struct TimeStepControlGrowthDampingFactor { static constexpr double value = 3.2; };
66struct TimeStepControlFileName { static constexpr auto value = "timesteps"; };
67struct MinTimeStepBeforeShuttingProblematicWellsInDays { static constexpr double value = 0.01; };
68struct MinTimeStepBasedOnNewtonIterations { static constexpr double value = 0.0; };
69
70} // namespace Opm::Parameters
71
72namespace Opm {
73
74struct StepReport;
75
76namespace detail {
77
78void logTimer(const AdaptiveSimulatorTimer& substepTimer);
79
80std::set<std::string> consistentlyFailingWells(const std::vector<StepReport>& sr);
81
83
84}
85
86 // AdaptiveTimeStepping
87 //---------------------
88 template<class TypeTag>
90 {
92 template <class Solver>
93 class SolutionTimeErrorSolverWrapper : public RelativeChangeInterface
94 {
95 const Solver& solver_;
96 public:
97 SolutionTimeErrorSolverWrapper(const Solver& solver)
98 : solver_(solver)
99 {}
100
102 double relativeChange() const
103 { return solver_.model().relativeChange(); }
104 };
105
106 template<class E>
107 void logException_(const E& exception, bool verbose)
108 {
109 if (verbose) {
110 std::string message;
111 message = "Caught Exception: ";
112 message += exception.what();
113 OpmLog::debug(message);
114 }
115 }
116
117 public:
119
121 AdaptiveTimeStepping(const UnitSystem& unitSystem,
122 const double max_next_tstep = -1.0,
123 const bool terminalOutput = true)
125 , restartFactor_(Parameters::Get<Parameters::SolverRestartFactor<Scalar>>()) // 0.33
126 , growthFactor_(Parameters::Get<Parameters::SolverGrowthFactor<Scalar>>()) // 2.0
127 , maxGrowth_(Parameters::Get<Parameters::SolverMaxGrowth<Scalar>>()) // 3.0
128 , maxTimeStep_(Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60) // 365.25
129 , minTimeStep_(unitSystem.to_si(UnitSystem::measure::time, Parameters::Get<Parameters::SolverMinTimeStep<Scalar>>())) // 1e-12;
130 , ignoreConvergenceFailure_(Parameters::Get<Parameters::SolverContinueOnConvergenceFailure>()) // false;
131 , solverRestartMax_(Parameters::Get<Parameters::SolverMaxRestarts>()) // 10
132 , solverVerbose_(Parameters::Get<Parameters::SolverVerbosity>() > 0 && terminalOutput) // 2
133 , timestepVerbose_(Parameters::Get<Parameters::TimeStepVerbosity>() > 0 && terminalOutput) // 2
134 , suggestedNextTimestep_((max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>() : max_next_tstep) * 24 * 60 * 60) // 1.0
135 , fullTimestepInitially_(Parameters::Get<Parameters::FullTimeStepInitially>()) // false
136 , timestepAfterEvent_(Parameters::Get<Parameters::TimeStepAfterEventInDays<Scalar>>() * 24 * 60 * 60) // 1e30
137 , useNewtonIteration_(false)
138 , minTimeStepBeforeShuttingProblematicWells_(Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day)
139
140 {
141 init_(unitSystem);
142 }
143
144
145
149 AdaptiveTimeStepping(double max_next_tstep,
150 const Tuning& tuning,
151 const UnitSystem& unitSystem,
152 const bool terminalOutput = true)
154 , restartFactor_(tuning.TSFCNV)
155 , growthFactor_(tuning.TFDIFF)
156 , maxGrowth_(tuning.TSFMAX)
157 , maxTimeStep_(tuning.TSMAXZ) // 365.25
158 , minTimeStep_(tuning.TSFMIN) // 0.1;
160 , solverRestartMax_(Parameters::Get<Parameters::SolverMaxRestarts>()) // 10
161 , solverVerbose_(Parameters::Get<Parameters::SolverVerbosity>() > 0 && terminalOutput) // 2
162 , timestepVerbose_(Parameters::Get<Parameters::TimeStepVerbosity>() > 0 && terminalOutput) // 2
163 , suggestedNextTimestep_(max_next_tstep <= 0 ? Parameters::Get<Parameters::InitialTimeStepInDays>() * 24 * 60 * 60 : max_next_tstep) // 1.0
164 , fullTimestepInitially_(Parameters::Get<Parameters::FullTimeStepInitially>()) // false
165 , timestepAfterEvent_(tuning.TMAXWC) // 1e30
166 , useNewtonIteration_(false)
167 , minTimeStepBeforeShuttingProblematicWells_(Parameters::Get<Parameters::MinTimeStepBeforeShuttingProblematicWellsInDays>() * unit::day)
168 {
169 init_(unitSystem);
170 }
171
172 static void registerParameters()
173 {
174 registerEclTimeSteppingParameters<Scalar>();
176 }
177
183 template <class Solver>
184 SimulatorReport step(const SimulatorTimer& simulatorTimer,
185 Solver& solver,
186 const bool isEvent,
187 const std::vector<int>* fipnum = nullptr,
188 const std::function<bool()> tuningUpdater = [](){return false;})
189 {
190 // Maybe update tuning
191 tuningUpdater();
192 SimulatorReport report;
193 const double timestep = simulatorTimer.currentStepLength();
194
195 // init last time step as a fraction of the given time step
196 if (suggestedNextTimestep_ < 0) {
198 }
199
201 suggestedNextTimestep_ = timestep;
202 }
203
204 // use seperate time step after event
205 if (isEvent && timestepAfterEvent_ > 0) {
207 }
208
209 auto& simulator = solver.model().simulator();
210 auto& problem = simulator.problem();
211
212 // create adaptive step timer with previously used sub step size
213 AdaptiveSimulatorTimer substepTimer(simulatorTimer, suggestedNextTimestep_, maxTimeStep_);
214
215 // counter for solver restarts
216 int restarts = 0;
217
218 // sub step time loop
219 while (!substepTimer.done()) {
220 // Maybe update tuning
221 // get current delta t
222 auto oldValue = suggestedNextTimestep_;
223 if (tuningUpdater()) {
224 // Use provideTimeStepEstimate to make we sure don't simulate longer than the report step is.
225 substepTimer.provideTimeStepEstimate(suggestedNextTimestep_);
226 suggestedNextTimestep_ = oldValue;
227 }
228 const double dt = substepTimer.currentStepLength();
229 if (timestepVerbose_) {
230 detail::logTimer(substepTimer);
231 }
232
233 SimulatorReportSingle substepReport;
234 std::string causeOfFailure;
235 try {
236 substepReport = solver.step(substepTimer);
237
238 if (solverVerbose_) {
239 // report number of linear iterations
240 OpmLog::debug("Overall linear iterations used: " + std::to_string(substepReport.total_linear_iterations));
241 }
242 }
243 catch (const TooManyIterations& e) {
244 substepReport = solver.failureReport();
245 causeOfFailure = "Solver convergence failure - Iteration limit reached";
246
247 logException_(e, solverVerbose_);
248 // since linearIterations is < 0 this will restart the solver
249 }
250 catch (const LinearSolverProblem& e) {
251 substepReport = solver.failureReport();
252 causeOfFailure = "Linear solver convergence failure";
253
254 logException_(e, solverVerbose_);
255 // since linearIterations is < 0 this will restart the solver
256 }
257 catch (const NumericalProblem& e) {
258 substepReport = solver.failureReport();
259 causeOfFailure = "Solver convergence failure - Numerical problem encountered";
260
261 logException_(e, solverVerbose_);
262 // since linearIterations is < 0 this will restart the solver
263 }
264 catch (const std::runtime_error& e) {
265 substepReport = solver.failureReport();
266
267 logException_(e, solverVerbose_);
268 // also catch linear solver not converged
269 }
270 catch (const Dune::ISTLError& e) {
271 substepReport = solver.failureReport();
272
273 logException_(e, solverVerbose_);
274 // also catch errors in ISTL AMG that occur when time step is too large
275 }
276 catch (const Dune::MatrixBlockError& e) {
277 substepReport = solver.failureReport();
278
279 logException_(e, solverVerbose_);
280 // this can be thrown by ISTL's ILU0 in block mode, yet is not an ISTLError
281 }
282
283 //Pass substep to eclwriter for summary output
284 simulator.problem().setSubStepReport(substepReport);
285
286 report += substepReport;
287
288 bool continue_on_uncoverged_solution = ignoreConvergenceFailure_ &&
289 !substepReport.converged &&
290 dt <= minTimeStep_;
291
292 if (continue_on_uncoverged_solution && solverVerbose_) {
293 const auto msg = fmt::format(
294 "Solver failed to converge but timestep "
295 "{} is smaller or equal to {}\n"
296 "which is the minimum threshold given "
297 "by option --solver-min-time-step\n",
298 dt, minTimeStep_
299 );
300 OpmLog::problem(msg);
301 }
302
303 if (substepReport.converged || continue_on_uncoverged_solution) {
304
305 // advance by current dt
306 ++substepTimer;
307
308 // create object to compute the time error, simply forwards the call to the model
309 SolutionTimeErrorSolverWrapper<Solver> relativeChange(solver);
310
311 // compute new time step estimate
312 const int iterations = useNewtonIteration_ ? substepReport.total_newton_iterations
313 : substepReport.total_linear_iterations;
314 double dtEstimate = timeStepControl_->computeTimeStepSize(dt, iterations, relativeChange,
315 substepTimer.simulationTimeElapsed());
316
317 assert(dtEstimate > 0);
318 // limit the growth of the timestep size by the growth factor
319 dtEstimate = std::min(dtEstimate, double(maxGrowth_ * dt));
320 assert(dtEstimate > 0);
321 // further restrict time step size growth after convergence problems
322 if (restarts > 0) {
323 dtEstimate = std::min(growthFactor_ * dt, dtEstimate);
324 // solver converged, reset restarts counter
325 restarts = 0;
326 }
327
328 // Further restrict time step size if we are in
329 // prediction mode with THP constraints.
330 if (solver.model().wellModel().hasTHPConstraints()) {
331 const double maxPredictionTHPTimestep = 16.0 * unit::day;
332 dtEstimate = std::min(dtEstimate, maxPredictionTHPTimestep);
333 }
334 assert(dtEstimate > 0);
335 if (timestepVerbose_) {
336 std::ostringstream ss;
337 substepReport.reportStep(ss);
338 OpmLog::info(ss.str());
339 }
340
341 // write data if outputWriter was provided
342 // if the time step is done we do not need
343 // to write it as this will be done by the simulator
344 // anyway.
345 if (!substepTimer.done()) {
346 if (fipnum) {
347 solver.computeFluidInPlace(*fipnum);
348 }
349 time::StopWatch perfTimer;
350 perfTimer.start();
351
352 problem.writeOutput();
353
354 report.success.output_write_time += perfTimer.secsSinceStart();
355 }
356
357 // set new time step length
358 substepTimer.provideTimeStepEstimate(dtEstimate);
359
360 report.success.converged = substepTimer.done();
361 substepTimer.setLastStepFailed(false);
362
363 }
364 else { // in case of no convergence
365 substepTimer.setLastStepFailed(true);
366
367 // If we have restarted (i.e. cut the timestep) too
368 // many times, we have failed and throw an exception.
369 if (restarts >= solverRestartMax_) {
370 const auto msg = fmt::format(
371 "Solver failed to converge after cutting timestep {} times.",
372 restarts
373 );
374 if (solverVerbose_) {
375 OpmLog::error(msg);
376 }
377 // Use throw directly to prevent file and line
378 throw TimeSteppingBreakdown{msg};
379 }
380
381 // The new, chopped timestep.
382 const double newTimeStep = restartFactor_ * dt;
383
384
385 // If we have restarted (i.e. cut the timestep) too
386 // much, we have failed and throw an exception.
387 if (newTimeStep < minTimeStep_) {
388 const auto msg = fmt::format(
389 "Solver failed to converge after cutting timestep to {}\n"
390 "which is the minimum threshold given by option --solver-min-time-step\n",
392 );
393 if (solverVerbose_) {
394 OpmLog::error(msg);
395 }
396 // Use throw directly to prevent file and line
397 throw TimeSteppingBreakdown{msg};
398 }
399
400 // Define utility function for chopping timestep.
401 auto chopTimestep = [&]() {
402 substepTimer.provideTimeStepEstimate(newTimeStep);
403 if (solverVerbose_) {
404 const auto msg = fmt::format(
405 "{}\nTimestep chopped to {} days\n",
406 causeOfFailure,
407 std::to_string(unit::convert::to(substepTimer.currentStepLength(), unit::day))
408 );
409 OpmLog::problem(msg);
410 }
411 ++restarts;
412 };
413
414 const double minimumChoppedTimestep = minTimeStepBeforeShuttingProblematicWells_;
415 if (newTimeStep > minimumChoppedTimestep) {
416 chopTimestep();
417 } else {
418 // We are below the threshold, and will check if there are any
419 // wells we should close rather than chopping again.
420 std::set<std::string> failing_wells = detail::consistentlyFailingWells(solver.model().stepReports());
421 if (failing_wells.empty()) {
422 // Found no wells to close, chop the timestep as above.
423 chopTimestep();
424 } else {
425 // Close all consistently failing wells.
426 int num_shut_wells = 0;
427 for (const auto& well : failing_wells) {
428 bool was_shut = solver.model().wellModel().forceShutWellByName(well, substepTimer.simulationTimeElapsed());
429 if (was_shut) {
430 ++num_shut_wells;
431 }
432 }
433 if (num_shut_wells == 0) {
434 // None of the problematic wells were shut.
435 // We must fall back to chopping again.
436 chopTimestep();
437 } else {
438 substepTimer.provideTimeStepEstimate(dt);
439 if (solverVerbose_) {
440 std::string msg;
441 msg = "\nProblematic well(s) were shut: ";
442 for (const auto& well : failing_wells) {
443 msg += well;
444 msg += " ";
445 }
446 msg += "(retrying timestep)\n";
447 OpmLog::problem(msg);
448 }
449 }
450 }
451 }
452 }
453 problem.setNextTimeStepSize(substepTimer.currentStepLength());
454 }
455
456 // store estimated time step for next reportStep
457 suggestedNextTimestep_ = substepTimer.currentStepLength();
458 if (timestepVerbose_) {
459 std::ostringstream ss;
460 substepTimer.report(ss);
461 ss << "Suggested next step size = " << unit::convert::to(suggestedNextTimestep_, unit::day) << " (days)" << std::endl;
462 OpmLog::debug(ss.str());
463 }
464
465 if (! std::isfinite(suggestedNextTimestep_)) { // check for NaN
466 suggestedNextTimestep_ = timestep;
467 }
468 return report;
469 }
470
474 double suggestedNextStep() const
475 { return suggestedNextTimestep_; }
476
477 void setSuggestedNextStep(const double x)
479
480 void updateTUNING(double max_next_tstep, const Tuning& tuning)
481 {
482 restartFactor_ = tuning.TSFCNV;
483 growthFactor_ = tuning.TFDIFF;
484 maxGrowth_ = tuning.TSFMAX;
485 maxTimeStep_ = tuning.TSMAXZ;
486 updateNEXTSTEP(max_next_tstep);
487 timestepAfterEvent_ = tuning.TMAXWC;
488 }
489
490 void updateNEXTSTEP(double max_next_tstep)
491 {
492 // \Note Only update next suggested step if TSINIT was explicitly set in TUNING or NEXTSTEP is active.
493 if (max_next_tstep > 0) {
494 suggestedNextTimestep_ = max_next_tstep;
495 }
496 }
497
498 template<class Serializer>
499 void serializeOp(Serializer& serializer)
500 {
501 serializer(timeStepControlType_);
502 switch (timeStepControlType_) {
504 allocAndSerialize<HardcodedTimeStepControl>(serializer);
505 break;
507 allocAndSerialize<PIDAndIterationCountTimeStepControl>(serializer);
508 break;
510 allocAndSerialize<SimpleIterationCountTimeStepControl>(serializer);
511 break;
513 allocAndSerialize<PIDTimeStepControl>(serializer);
514 break;
515 }
516 serializer(restartFactor_);
517 serializer(growthFactor_);
518 serializer(maxGrowth_);
519 serializer(maxTimeStep_);
520 serializer(minTimeStep_);
521 serializer(ignoreConvergenceFailure_);
522 serializer(solverRestartMax_);
523 serializer(solverVerbose_);
524 serializer(timestepVerbose_);
525 serializer(suggestedNextTimestep_);
526 serializer(fullTimestepInitially_);
527 serializer(timestepAfterEvent_);
528 serializer(useNewtonIteration_);
530 }
531
533 {
534 return serializationTestObject_<HardcodedTimeStepControl>();
535 }
536
538 {
539 return serializationTestObject_<PIDTimeStepControl>();
540 }
541
543 {
544 return serializationTestObject_<PIDAndIterationCountTimeStepControl>();
545 }
546
548 {
549 return serializationTestObject_<SimpleIterationCountTimeStepControl>();
550 }
551
553 {
557 return false;
558 }
559
560 bool result = false;
561 switch (timeStepControlType_) {
563 result = castAndComp<HardcodedTimeStepControl>(rhs);
564 break;
566 result = castAndComp<PIDAndIterationCountTimeStepControl>(rhs);
567 break;
569 result = castAndComp<SimpleIterationCountTimeStepControl>(rhs);
570 break;
572 result = castAndComp<PIDTimeStepControl>(rhs);
573 break;
574 }
575
576 return result &&
577 this->restartFactor_ == rhs.restartFactor_ &&
578 this->growthFactor_ == rhs.growthFactor_ &&
579 this->maxGrowth_ == rhs.maxGrowth_ &&
580 this->maxTimeStep_ == rhs.maxTimeStep_ &&
581 this->minTimeStep_ == rhs.minTimeStep_ &&
584 this->solverVerbose_ == rhs.solverVerbose_ &&
590 }
591
592 private:
593 template<class Controller>
594 static AdaptiveTimeStepping<TypeTag> serializationTestObject_()
595 {
597
598 result.restartFactor_ = 1.0;
599 result.growthFactor_ = 2.0;
600 result.maxGrowth_ = 3.0;
601 result.maxTimeStep_ = 4.0;
602 result.minTimeStep_ = 5.0;
603 result.ignoreConvergenceFailure_ = true;
604 result.solverRestartMax_ = 6;
605 result.solverVerbose_ = true;
606 result.timestepVerbose_ = true;
607 result.suggestedNextTimestep_ = 7.0;
608 result.fullTimestepInitially_ = true;
609 result.useNewtonIteration_ = true;
611 result.timeStepControlType_ = Controller::Type;
612 result.timeStepControl_ = std::make_unique<Controller>(Controller::serializationTestObject());
613
614 return result;
615 }
616 template<class T, class Serializer>
617 void allocAndSerialize(Serializer& serializer)
618 {
619 if (!serializer.isSerializing()) {
620 timeStepControl_ = std::make_unique<T>();
621 }
622 serializer(*static_cast<T*>(timeStepControl_.get()));
623 }
624
625 template<class T>
626 bool castAndComp(const AdaptiveTimeStepping<TypeTag>& Rhs) const
627 {
628 const T* lhs = static_cast<const T*>(timeStepControl_.get());
629 const T* rhs = static_cast<const T*>(Rhs.timeStepControl_.get());
630 return *lhs == *rhs;
631 }
632
633 protected:
634 void init_(const UnitSystem& unitSystem)
635 {
636 // valid are "pid" and "pid+iteration"
637 std::string control = Parameters::Get<Parameters::TimeStepControl>(); // "pid"
638
639 const double tol = Parameters::Get<Parameters::TimeStepControlTolerance>(); // 1e-1
640 if (control == "pid") {
641 timeStepControl_ = std::make_unique<PIDTimeStepControl>(tol);
643 }
644 else if (control == "pid+iteration") {
645 const int iterations = Parameters::Get<Parameters::TimeStepControlTargetIterations>(); // 30
646 const double decayDampingFactor = Parameters::Get<Parameters::TimeStepControlDecayDampingFactor>(); // 1.0
647 const double growthDampingFactor = Parameters::Get<Parameters::TimeStepControlGrowthDampingFactor>(); // 3.2
648 timeStepControl_ = std::make_unique<PIDAndIterationCountTimeStepControl>(iterations, decayDampingFactor, growthDampingFactor, tol);
650 }
651 else if (control == "pid+newtoniteration") {
652 const int iterations = Parameters::Get<Parameters::TimeStepControlTargetNewtonIterations>(); // 8
653 const double decayDampingFactor = Parameters::Get<Parameters::TimeStepControlDecayDampingFactor>(); // 1.0
654 const double growthDampingFactor = Parameters::Get<Parameters::TimeStepControlGrowthDampingFactor>(); // 3.2
655 const double nonDimensionalMinTimeStepIterations = Parameters::Get<Parameters::MinTimeStepBasedOnNewtonIterations>(); // 0.0 by default
656 // the min time step can be reduced by the newton iteration numbers
657 double minTimeStepReducedByIterations = unitSystem.to_si(UnitSystem::measure::time, nonDimensionalMinTimeStepIterations);
658 timeStepControl_ = std::make_unique<PIDAndIterationCountTimeStepControl>(iterations, decayDampingFactor,
659 growthDampingFactor, tol, minTimeStepReducedByIterations);
661 useNewtonIteration_ = true;
662 }
663 else if (control == "iterationcount") {
664 const int iterations = Parameters::Get<Parameters::TimeStepControlTargetIterations>(); // 30
665 const double decayrate = Parameters::Get<Parameters::TimeStepControlDecayRate>(); // 0.75
666 const double growthrate = Parameters::Get<Parameters::TimeStepControlGrowthRate>(); // 1.25
667 timeStepControl_ = std::make_unique<SimpleIterationCountTimeStepControl>(iterations, decayrate, growthrate);
669 }
670 else if (control == "newtoniterationcount") {
671 const int iterations = Parameters::Get<Parameters::TimeStepControlTargetNewtonIterations>(); // 8
672 const double decayrate = Parameters::Get<Parameters::TimeStepControlDecayRate>(); // 0.75
673 const double growthrate = Parameters::Get<Parameters::TimeStepControlGrowthRate>(); // 1.25
674 timeStepControl_ = std::make_unique<SimpleIterationCountTimeStepControl>(iterations, decayrate, growthrate);
675 useNewtonIteration_ = true;
677 }
678 else if (control == "hardcoded") {
679 const std::string filename = Parameters::Get<Parameters::TimeStepControlFileName>(); // "timesteps"
680 timeStepControl_ = std::make_unique<HardcodedTimeStepControl>(filename);
682 }
683 else
684 OPM_THROW(std::runtime_error,
685 "Unsupported time step control selected " + control);
686
687 // make sure growth factor is something reasonable
688 assert(growthFactor_ >= 1.0);
689 }
690
691 using TimeStepController = std::unique_ptr<TimeStepControlInterface>;
692
697 double maxGrowth_;
709 };
710}
711
712#endif // OPM_ADAPTIVE_TIME_STEPPING_HPP
Defines a type tags and some fundamental properties all models.
Simulation timer for adaptive time stepping.
Definition: AdaptiveSimulatorTimer.hpp:41
Definition: AdaptiveTimeStepping.hpp:90
double minTimeStep_
minimal allowed time step size before throwing
Definition: AdaptiveTimeStepping.hpp:699
double suggestedNextTimestep_
suggested size of next timestep
Definition: AdaptiveTimeStepping.hpp:704
bool solverVerbose_
solver verbosity
Definition: AdaptiveTimeStepping.hpp:702
void init_(const UnitSystem &unitSystem)
Definition: AdaptiveTimeStepping.hpp:634
void setSuggestedNextStep(const double x)
Definition: AdaptiveTimeStepping.hpp:477
double suggestedNextStep() const
Returns the simulator report for the failed substeps of the last report step.
Definition: AdaptiveTimeStepping.hpp:474
double growthFactor_
factor to multiply time step when solver recovered from failed convergence
Definition: AdaptiveTimeStepping.hpp:696
AdaptiveTimeStepping(const UnitSystem &unitSystem, const double max_next_tstep=-1.0, const bool terminalOutput=true)
contructor taking parameter object
Definition: AdaptiveTimeStepping.hpp:121
void serializeOp(Serializer &serializer)
Definition: AdaptiveTimeStepping.hpp:499
void updateTUNING(double max_next_tstep, const Tuning &tuning)
Definition: AdaptiveTimeStepping.hpp:480
static void registerParameters()
Definition: AdaptiveTimeStepping.hpp:172
TimeStepController timeStepControl_
time step control object
Definition: AdaptiveTimeStepping.hpp:694
static AdaptiveTimeStepping< TypeTag > serializationTestObjectSimple()
Definition: AdaptiveTimeStepping.hpp:547
bool fullTimestepInitially_
beginning with the size of the time step from data file
Definition: AdaptiveTimeStepping.hpp:705
AdaptiveTimeStepping(double max_next_tstep, const Tuning &tuning, const UnitSystem &unitSystem, const bool terminalOutput=true)
contructor taking parameter object
Definition: AdaptiveTimeStepping.hpp:149
std::unique_ptr< TimeStepControlInterface > TimeStepController
Definition: AdaptiveTimeStepping.hpp:691
double timestepAfterEvent_
suggested size of timestep after an event
Definition: AdaptiveTimeStepping.hpp:706
bool useNewtonIteration_
use newton iteration count for adaptive time step control
Definition: AdaptiveTimeStepping.hpp:707
bool timestepVerbose_
timestep verbosity
Definition: AdaptiveTimeStepping.hpp:703
int solverRestartMax_
how many restart of solver are allowed
Definition: AdaptiveTimeStepping.hpp:701
TimeStepControlType timeStepControlType_
type of time step control object
Definition: AdaptiveTimeStepping.hpp:693
double minTimeStepBeforeShuttingProblematicWells_
Definition: AdaptiveTimeStepping.hpp:708
bool operator==(const AdaptiveTimeStepping< TypeTag > &rhs)
Definition: AdaptiveTimeStepping.hpp:552
SimulatorReport step(const SimulatorTimer &simulatorTimer, Solver &solver, const bool isEvent, const std::vector< int > *fipnum=nullptr, const std::function< bool()> tuningUpdater=[](){return false;})
step method that acts like the solver::step method in a sub cycle of time steps
Definition: AdaptiveTimeStepping.hpp:184
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPIDIt()
Definition: AdaptiveTimeStepping.hpp:542
void updateNEXTSTEP(double max_next_tstep)
Definition: AdaptiveTimeStepping.hpp:490
bool ignoreConvergenceFailure_
continue instead of stop when minimum time step is reached
Definition: AdaptiveTimeStepping.hpp:700
double maxTimeStep_
maximal allowed time step size in days
Definition: AdaptiveTimeStepping.hpp:698
static AdaptiveTimeStepping< TypeTag > serializationTestObjectHardcoded()
Definition: AdaptiveTimeStepping.hpp:532
double restartFactor_
factor to multiply time step with when solver fails to converge
Definition: AdaptiveTimeStepping.hpp:695
double maxGrowth_
factor that limits the maximum growth of a time step
Definition: AdaptiveTimeStepping.hpp:697
static AdaptiveTimeStepping< TypeTag > serializationTestObjectPID()
Definition: AdaptiveTimeStepping.hpp:537
Definition: TimeStepControlInterface.hpp:32
Definition: SimulatorTimer.hpp:39
double currentStepLength() const override
Definition: blackoilnewtonmethodparams.hpp:31
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:185
void logTimer(const AdaptiveSimulatorTimer &substepTimer)
void registerAdaptiveParameters()
std::set< std::string > consistentlyFailingWells(const std::vector< StepReport > &sr)
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
TimeStepControlType
Definition: TimeStepControl.hpp:31
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: AdaptiveTimeStepping.hpp:57
static constexpr bool value
Definition: AdaptiveTimeStepping.hpp:57
Definition: AdaptiveTimeStepping.hpp:56
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:56
Definition: AdaptiveTimeStepping.hpp:68
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:68
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:67
Definition: AdaptiveTimeStepping.hpp:52
static constexpr bool value
Definition: AdaptiveTimeStepping.hpp:52
Definition: AdaptiveTimeStepping.hpp:53
static constexpr int value
Definition: AdaptiveTimeStepping.hpp:53
Definition: AdaptiveTimeStepping.hpp:54
static constexpr int value
Definition: AdaptiveTimeStepping.hpp:54
Definition: AdaptiveTimeStepping.hpp:64
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:64
Definition: AdaptiveTimeStepping.hpp:62
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:62
Definition: AdaptiveTimeStepping.hpp:66
static constexpr auto value
Definition: AdaptiveTimeStepping.hpp:66
Definition: AdaptiveTimeStepping.hpp:65
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:65
Definition: AdaptiveTimeStepping.hpp:63
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:63
Definition: AdaptiveTimeStepping.hpp:58
static constexpr auto value
Definition: AdaptiveTimeStepping.hpp:58
Definition: AdaptiveTimeStepping.hpp:60
static constexpr int value
Definition: AdaptiveTimeStepping.hpp:60
Definition: AdaptiveTimeStepping.hpp:61
static constexpr int value
Definition: AdaptiveTimeStepping.hpp:61
Definition: AdaptiveTimeStepping.hpp:59
static constexpr double value
Definition: AdaptiveTimeStepping.hpp:59
Definition: AdaptiveTimeStepping.hpp:55
static constexpr int value
Definition: AdaptiveTimeStepping.hpp:55
Definition: SimulatorReport.hpp:100
SimulatorReportSingle success
Definition: SimulatorReport.hpp:101
bool converged
Definition: SimulatorReport.hpp:54
double output_write_time
Definition: SimulatorReport.hpp:45
Definition: ConvergenceReport.hpp:338