3 #ifndef OPM_ADAPTIVE_TIME_STEPPING_EBOS_HPP
4 #define OPM_ADAPTIVE_TIME_STEPPING_EBOS_HPP
9 #include <opm/simulators/timestepping/SimulatorReport.hpp>
10 #include <opm/grid/utility/StopWatch.hpp>
11 #include <opm/common/OpmLog/OpmLog.hpp>
12 #include <opm/common/utility/parameters/ParameterGroup.hpp>
13 #include <opm/common/ErrorMacros.hpp>
14 #include <opm/simulators/timestepping/SimulatorTimer.hpp>
15 #include <opm/simulators/timestepping/AdaptiveSimulatorTimer.hpp>
16 #include <opm/simulators/timestepping/TimeStepControlInterface.hpp>
17 #include <opm/simulators/timestepping/TimeStepControl.hpp>
18 #include <opm/core/props/phaseUsageFromDeck.hpp>
19 #include <opm/common/Exceptions.hpp>
21 namespace Opm::Properties {
27 template<
class TypeTag,
class MyTypeTag>
29 using type = UndefinedProperty;
31 template<
class TypeTag,
class MyTypeTag>
33 using type = UndefinedProperty;
35 template<
class TypeTag,
class MyTypeTag>
37 using type = UndefinedProperty;
39 template<
class TypeTag,
class MyTypeTag>
41 using type = UndefinedProperty;
43 template<
class TypeTag,
class MyTypeTag>
45 using type = UndefinedProperty;
47 template<
class TypeTag,
class MyTypeTag>
49 using type = UndefinedProperty;
51 template<
class TypeTag,
class MyTypeTag>
53 using type = UndefinedProperty;
55 template<
class TypeTag,
class MyTypeTag>
57 using type = UndefinedProperty;
59 template<
class TypeTag,
class MyTypeTag>
61 using type = UndefinedProperty;
63 template<
class TypeTag,
class MyTypeTag>
65 using type = UndefinedProperty;
67 template<
class TypeTag,
class MyTypeTag>
69 using type = UndefinedProperty;
71 template<
class TypeTag,
class MyTypeTag>
73 using type = UndefinedProperty;
75 template<
class TypeTag,
class MyTypeTag>
77 using type = UndefinedProperty;
79 template<
class TypeTag,
class MyTypeTag>
81 using type = UndefinedProperty;
83 template<
class TypeTag,
class MyTypeTag>
85 using type = UndefinedProperty;
87 template<
class TypeTag,
class MyTypeTag>
89 using type = UndefinedProperty;
91 template<
class TypeTag,
class MyTypeTag>
93 using type = UndefinedProperty;
95 template<
class TypeTag,
class MyTypeTag>
97 using type = UndefinedProperty;
99 template<
class TypeTag,
class MyTypeTag>
101 using type = UndefinedProperty;
103 template<
class TypeTag,
class MyTypeTag>
105 using type = UndefinedProperty;
107 template<
class TypeTag,
class MyTypeTag>
109 using type = UndefinedProperty;
111 template<
class TypeTag,
class MyTypeTag>
113 using type = UndefinedProperty;
115 template<
class TypeTag,
class MyTypeTag>
117 using type = UndefinedProperty;
120 template<
class TypeTag>
122 using type = GetPropType<TypeTag, Scalar>;
123 static constexpr type value = 0.33;
125 template<
class TypeTag>
127 using type = GetPropType<TypeTag, Scalar>;
128 static constexpr type value = 2.0;
130 template<
class TypeTag>
132 using type = GetPropType<TypeTag, Scalar>;
133 static constexpr type value = 3.0;
135 template<
class TypeTag>
137 using type = GetPropType<TypeTag, Scalar>;
138 static constexpr type value = 365.0;
140 template<
class TypeTag>
142 using type = GetPropType<TypeTag, Scalar>;
143 static constexpr type value = 1.0e-12;
145 template<
class TypeTag>
147 static constexpr
bool value =
false;
149 template<
class TypeTag>
151 static constexpr
int value = 10;
153 template<
class TypeTag>
155 static constexpr
int value = 1;
157 template<
class TypeTag>
159 static constexpr
int value = 1;
161 template<
class TypeTag>
163 using type = GetPropType<TypeTag, Scalar>;
164 static constexpr type value = 1.0;
166 template<
class TypeTag>
168 static constexpr
bool value =
false;
170 template<
class TypeTag>
172 using type = GetPropType<TypeTag, Scalar>;
173 static constexpr type value = -1.0;
175 template<
class TypeTag>
177 static constexpr
auto value =
"pid+newtoniteration";
179 template<
class TypeTag>
181 using type = GetPropType<TypeTag, Scalar>;
182 static constexpr type value = 1e-1;
184 template<
class TypeTag>
186 static constexpr
int value = 30;
188 template<
class TypeTag>
190 static constexpr
int value = 8;
192 template<
class TypeTag>
194 using type = GetPropType<TypeTag, Scalar>;
195 static constexpr type value = 0.75;
197 template<
class TypeTag>
199 using type = GetPropType<TypeTag, Scalar>;
200 static constexpr type value = 1.25;
202 template<
class TypeTag>
204 using type = GetPropType<TypeTag, Scalar>;
205 static constexpr type value = 1.0;
207 template<
class TypeTag>
209 using type = GetPropType<TypeTag, Scalar>;
210 static constexpr type value = 3.2;
212 template<
class TypeTag>
214 static constexpr
auto value =
"timesteps";
216 template<
class TypeTag>
218 using type = GetPropType<TypeTag, Scalar>;
219 static constexpr type value = 0.01;
222 template<
class TypeTag>
224 using type = GetPropType<TypeTag, Scalar>;
225 static constexpr type value = 0.0;
235 template<
class TypeTag>
238 template <
class Solver>
241 const Solver& solver_;
243 SolutionTimeErrorSolverWrapperEbos(
const Solver& solver)
248 double relativeChange()
const
249 {
return solver_.model().relativeChange(); }
253 void logException_(
const E& exception,
bool verbose)
257 message =
"Caught Exception: ";
258 message += exception.what();
259 OpmLog::debug(message);
266 const bool terminalOutput =
true)
268 ,
restartFactor_(EWOMS_GET_PARAM(TypeTag, double, SolverRestartFactor))
269 ,
growthFactor_(EWOMS_GET_PARAM(TypeTag, double, SolverGrowthFactor))
270 ,
maxGrowth_(EWOMS_GET_PARAM(TypeTag, double, SolverMaxGrowth))
271 ,
maxTimeStep_(EWOMS_GET_PARAM(TypeTag, double, SolverMaxTimeStepInDays)*24*60*60)
272 ,
minTimeStep_(unitSystem.to_si(UnitSystem::measure::time, EWOMS_GET_PARAM(TypeTag, double, SolverMinTimeStep)))
275 ,
solverVerbose_(EWOMS_GET_PARAM(TypeTag, int, SolverVerbosity) > 0 && terminalOutput)
276 ,
timestepVerbose_(EWOMS_GET_PARAM(TypeTag, int, TimeStepVerbosity) > 0 && terminalOutput)
281 , minTimeStepBeforeShuttingProblematicWells_(EWOMS_GET_PARAM(TypeTag, double, MinTimeStepBeforeShuttingProblematicWellsInDays)*unit::day)
293 const UnitSystem& unitSystem,
294 const bool terminalOutput =
true)
303 ,
solverVerbose_(EWOMS_GET_PARAM(TypeTag, int, SolverVerbosity) > 0 && terminalOutput)
304 ,
timestepVerbose_(EWOMS_GET_PARAM(TypeTag, int, TimeStepVerbosity) > 0 && terminalOutput)
309 , minTimeStepBeforeShuttingProblematicWells_(EWOMS_GET_PARAM(TypeTag, double, MinTimeStepBeforeShuttingProblematicWellsInDays)*unit::day)
314 static void registerParameters()
317 EWOMS_REGISTER_PARAM(TypeTag,
double, SolverRestartFactor,
318 "The factor time steps are elongated after restarts");
319 EWOMS_REGISTER_PARAM(TypeTag,
double, SolverGrowthFactor,
320 "The factor time steps are elongated after a successful substep");
321 EWOMS_REGISTER_PARAM(TypeTag,
double, SolverMaxGrowth,
322 "The maximum factor time steps are elongated after a report step");
323 EWOMS_REGISTER_PARAM(TypeTag,
double, SolverMaxTimeStepInDays,
324 "The maximum size of a time step in days");
325 EWOMS_REGISTER_PARAM(TypeTag,
double, SolverMinTimeStep,
326 "The minimum size of a time step in days for field and metric and hours for lab. If a step cannot converge without getting cut below this step size the simulator will stop");
327 EWOMS_REGISTER_PARAM(TypeTag,
bool, SolverContinueOnConvergenceFailure,
328 "Continue instead of stop when minimum solver time step is reached");
329 EWOMS_REGISTER_PARAM(TypeTag,
int, SolverMaxRestarts,
330 "The maximum number of breakdowns before a substep is given up and the simulator is terminated");
331 EWOMS_REGISTER_PARAM(TypeTag,
int, SolverVerbosity,
332 "Specify the \"chattiness\" of the non-linear solver itself");
333 EWOMS_REGISTER_PARAM(TypeTag,
int, TimeStepVerbosity,
334 "Specify the \"chattiness\" during the time integration");
335 EWOMS_REGISTER_PARAM(TypeTag,
double, InitialTimeStepInDays,
336 "The size of the initial time step in days");
337 EWOMS_REGISTER_PARAM(TypeTag,
bool, FullTimeStepInitially,
338 "Always attempt to finish a report step using a single substep");
339 EWOMS_REGISTER_PARAM(TypeTag,
double, TimeStepAfterEventInDays,
340 "Time step size of the first time step after an event occurs during the simulation in days");
341 EWOMS_REGISTER_PARAM(TypeTag, std::string, TimeStepControl,
342 "The algorithm used to determine time-step sizes. valid options are: 'pid' (default), 'pid+iteration', 'pid+newtoniteration', 'iterationcount', 'newtoniterationcount' and 'hardcoded'");
343 EWOMS_REGISTER_PARAM(TypeTag,
double, TimeStepControlTolerance,
344 "The tolerance used by the time step size control algorithm");
345 EWOMS_REGISTER_PARAM(TypeTag,
int, TimeStepControlTargetIterations,
346 "The number of linear iterations which the time step control scheme should aim for (if applicable)");
347 EWOMS_REGISTER_PARAM(TypeTag,
int, TimeStepControlTargetNewtonIterations,
348 "The number of Newton iterations which the time step control scheme should aim for (if applicable)");
349 EWOMS_REGISTER_PARAM(TypeTag,
double, TimeStepControlDecayRate,
350 "The decay rate of the time step size of the number of target iterations is exceeded");
351 EWOMS_REGISTER_PARAM(TypeTag,
double, TimeStepControlGrowthRate,
352 "The growth rate of the time step size of the number of target iterations is undercut");
353 EWOMS_REGISTER_PARAM(TypeTag,
double, TimeStepControlDecayDampingFactor,
354 "The decay rate of the time step decrease when the target iterations is exceeded");
355 EWOMS_REGISTER_PARAM(TypeTag,
double, TimeStepControlGrowthDampingFactor,
356 "The growth rate of the time step increase when the target iterations is undercut");
357 EWOMS_REGISTER_PARAM(TypeTag, std::string, TimeStepControlFileName,
358 "The name of the file which contains the hardcoded time steps sizes");
359 EWOMS_REGISTER_PARAM(TypeTag,
double, MinTimeStepBeforeShuttingProblematicWellsInDays,
360 "The minimum time step size in days for which problematic wells are not shut");
361 EWOMS_REGISTER_PARAM(TypeTag,
double, MinTimeStepBasedOnNewtonIterations,
362 "The minimum time step size (in days for field and metric unit and hours for lab unit) can be reduced to based on newton iteration counts");
368 template <
class Solver>
372 const std::vector<int>* fipnum =
nullptr)
391 auto& ebosSimulator = solver.model().ebosSimulator();
392 auto& ebosProblem = ebosSimulator.problem();
401 while (!substepTimer.
done()) {
405 logTimer(substepTimer);
409 std::string causeOfFailure =
"";
411 substepReport = solver.step(substepTimer);
414 OpmLog::debug(
"Overall linear iterations used: " + std::to_string(substepReport.total_linear_iterations));
417 catch (
const TooManyIterations& e) {
418 substepReport = solver.failureReport();
419 causeOfFailure =
"Solver convergence failure - Iteration limit reached";
424 catch (
const LinearSolverProblem& e) {
425 substepReport = solver.failureReport();
426 causeOfFailure =
"Linear solver convergence failure";
431 catch (
const NumericalIssue& e) {
432 substepReport = solver.failureReport();
433 causeOfFailure =
"Solver convergence failure - Numerical problem encountered";
438 catch (
const std::runtime_error& e) {
439 substepReport = solver.failureReport();
444 catch (
const Dune::ISTLError& e) {
445 substepReport = solver.failureReport();
450 catch (
const Dune::MatrixBlockError& e) {
451 substepReport = solver.failureReport();
457 report += substepReport;
461 if (continue_on_uncoverged_solution) {
462 const auto msg = std::string(
"Solver failed to converge but timestep ")
463 + std::to_string(dt) +
" is smaller or equal to "
464 + std::to_string(
minTimeStep_) +
"\n which is the minimum threshold given"
465 +
"by option --solver-min-time-step= \n";
471 if (substepReport.converged || continue_on_uncoverged_solution) {
477 SolutionTimeErrorSolverWrapperEbos<Solver> relativeChange(solver);
481 : substepReport.total_linear_iterations;
482 double dtEstimate =
timeStepControl_->computeTimeStepSize(dt, iterations, relativeChange,
485 assert(dtEstimate > 0);
487 dtEstimate = std::min(dtEstimate,
double(
maxGrowth_ * dt));
488 assert(dtEstimate > 0);
498 if (solver.model().wellModel().hasTHPConstraints()) {
499 const double maxPredictionTHPTimestep = 16.0 * unit::day;
500 dtEstimate = std::min(dtEstimate, maxPredictionTHPTimestep);
502 assert(dtEstimate > 0);
504 std::ostringstream ss;
506 OpmLog::info(ss.str());
513 if (!substepTimer.
done()) {
515 solver.computeFluidInPlace(*fipnum);
517 time::StopWatch perfTimer;
520 ebosProblem.writeOutput();
522 report.success.output_write_time += perfTimer.secsSinceStart();
528 report.success.converged = substepTimer.
done();
538 const auto msg = std::string(
"Solver failed to converge after cutting timestep ")
539 + std::to_string(restarts) +
" times.";
543 OPM_THROW_NOLOG(NumericalIssue, msg);
553 const auto msg = std::string(
"Solver failed to converge after cutting timestep to ")
554 + std::to_string(
minTimeStep_) +
"\n which is the minimum threshold given"
555 +
"by option --solver-min-time-step= \n";
559 OPM_THROW_NOLOG(NumericalIssue, msg);
563 auto chopTimestep = [&]() {
567 msg = causeOfFailure +
"\nTimestep chopped to "
568 + std::to_string(unit::convert::to(substepTimer.
currentStepLength(), unit::day)) +
" days\n";
569 OpmLog::problem(msg);
574 const double minimumChoppedTimestep = minTimeStepBeforeShuttingProblematicWells_;
575 if (newTimeStep > minimumChoppedTimestep) {
580 std::set<std::string> failing_wells = consistentlyFailingWells(solver.model().stepReports());
581 if (failing_wells.empty()) {
586 int num_shut_wells = 0;
587 for (
const auto& well : failing_wells) {
588 bool was_shut = solver.model().wellModel().forceShutWellByNameIfPredictionMode(well, substepTimer.
simulationTimeElapsed());
593 if (num_shut_wells == 0) {
601 msg =
"\nProblematic well(s) were shut: ";
602 for (
const auto& well : failing_wells) {
606 msg +=
"(retrying timestep)\n";
607 OpmLog::problem(msg);
619 std::ostringstream ss;
621 ss <<
"Suggested next step size = " << unit::convert::to(
suggestedNextTimestep_, unit::day) <<
" (days)" << std::endl;
622 OpmLog::debug(ss.str());
637 void setSuggestedNextStep(
const double x)
640 void updateTUNING(
const Tuning& tuning)
652 void init_(
const UnitSystem& unitSystem)
655 std::string control = EWOMS_GET_PARAM(TypeTag, std::string, TimeStepControl);
657 const double tol = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlTolerance);
658 if (control ==
"pid") {
661 else if (control ==
"pid+iteration") {
662 const int iterations = EWOMS_GET_PARAM(TypeTag,
int, TimeStepControlTargetIterations);
663 const double decayDampingFactor = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlDecayDampingFactor);
664 const double growthDampingFactor = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlGrowthDampingFactor);
665 timeStepControl_ = TimeStepControlType(
new PIDAndIterationCountTimeStepControl(iterations, decayDampingFactor, growthDampingFactor, tol));
667 else if (control ==
"pid+newtoniteration") {
668 const int iterations = EWOMS_GET_PARAM(TypeTag,
int, TimeStepControlTargetNewtonIterations);
669 const double decayDampingFactor = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlDecayDampingFactor);
670 const double growthDampingFactor = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlGrowthDampingFactor);
671 const double nonDimensionalMinTimeStepIterations = EWOMS_GET_PARAM(TypeTag,
double, MinTimeStepBasedOnNewtonIterations);
673 double minTimeStepReducedByIterations = unitSystem.to_si(UnitSystem::measure::time, nonDimensionalMinTimeStepIterations);
674 timeStepControl_ = TimeStepControlType(
new PIDAndIterationCountTimeStepControl(iterations, decayDampingFactor,
675 growthDampingFactor, tol, minTimeStepReducedByIterations));
678 else if (control ==
"iterationcount") {
679 const int iterations = EWOMS_GET_PARAM(TypeTag,
int, TimeStepControlTargetIterations);
680 const double decayrate = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlDecayRate);
681 const double growthrate = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlGrowthRate);
682 timeStepControl_ = TimeStepControlType(
new SimpleIterationCountTimeStepControl(iterations, decayrate, growthrate));
684 else if (control ==
"newtoniterationcount") {
685 const int iterations = EWOMS_GET_PARAM(TypeTag,
int, TimeStepControlTargetNewtonIterations);
686 const double decayrate = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlDecayRate);
687 const double growthrate = EWOMS_GET_PARAM(TypeTag,
double, TimeStepControlGrowthRate);
688 timeStepControl_ = TimeStepControlType(
new SimpleIterationCountTimeStepControl(iterations, decayrate, growthrate));
691 else if (control ==
"hardcoded") {
692 const std::string filename = EWOMS_GET_PARAM(TypeTag, std::string, TimeStepControlFileName);
693 timeStepControl_ = TimeStepControlType(
new HardcodedTimeStepControl(filename));
697 OPM_THROW(std::runtime_error,
"Unsupported time step control selected "<< control);
703 template <
class ProblemType>
704 std::set<std::string> consistentlyFailingWells(
const std::vector<ProblemType>& sr)
708 std::ostringstream msg;
709 msg <<
" Excessive chopping detected in report step "
710 << sr.back().report_step <<
", substep " << sr.back().current_step <<
"\n";
712 std::set<std::string> failing_wells;
716 if(sr.back().report.empty())
717 return failing_wells;
719 const auto& wfs = sr.back().report.back().wellFailures();
720 for (
const auto& wf : wfs) {
721 msg <<
" Well that failed: " << wf.wellName() <<
"\n";
724 OpmLog::debug(msg.str());
727 const int num_steps = 3;
728 const int rep_step = sr.back().report_step;
729 const int sub_step = sr.back().current_step;
730 const int sr_size = sr.size();
731 if (sr_size >= num_steps) {
732 for (
const auto& wf : wfs) {
733 failing_wells.insert(wf.wellName());
735 for (
int s = 1; s < num_steps; ++s) {
736 const auto& srep = sr[sr_size - 1 - s];
739 if (srep.report_step != rep_step || srep.current_step != sub_step) {
743 std::set<std::string> failing_wells_step;
744 for (
const auto& wf : srep.report.back().wellFailures()) {
745 if (failing_wells.count(wf.wellName()) > 0) {
746 failing_wells_step.insert(wf.wellName());
749 failing_wells.swap(failing_wells_step);
752 return failing_wells;
755 typedef std::unique_ptr<TimeStepControlInterface> TimeStepControlType;
771 double minTimeStepBeforeShuttingProblematicWells_;
Simulation timer for adaptive time stepping.
Definition: AdaptiveSimulatorTimer.hpp:41
bool done() const
Definition: AdaptiveSimulatorTimer.cpp:126
double currentStepLength() const
Definition: AdaptiveSimulatorTimer.cpp:108
void provideTimeStepEstimate(const double dt_estimate)
provide and estimate for new time step size
Definition: AdaptiveSimulatorTimer.cpp:72
void report(std::ostream &os) const
report start and end time as well as used steps so far
Definition: AdaptiveSimulatorTimer.cpp:153
double simulationTimeElapsed() const
Definition: AdaptiveSimulatorTimer.cpp:124
void setLastStepFailed(bool lastStepFailed)
tell the timestepper whether timestep failed or not
Definition: AdaptiveSimulatorTimer.hpp:104
Definition: AdaptiveTimeSteppingEbos.hpp:237
double maxTimeStep_
maximal allowed time step size in days
Definition: AdaptiveTimeSteppingEbos.hpp:761
double growthFactor_
factor to multiply time step when solver recovered from failed convergence
Definition: AdaptiveTimeSteppingEbos.hpp:759
double minTimeStep_
minimal allowed time step size before throwing
Definition: AdaptiveTimeSteppingEbos.hpp:762
bool fullTimestepInitially_
beginning with the size of the time step from data file
Definition: AdaptiveTimeSteppingEbos.hpp:768
double suggestedNextStep() const
Returns the simulator report for the failed substeps of the last report step.
Definition: AdaptiveTimeSteppingEbos.hpp:634
int solverRestartMax_
how many restart of solver are allowed
Definition: AdaptiveTimeSteppingEbos.hpp:764
double timestepAfterEvent_
suggested size of timestep after an event
Definition: AdaptiveTimeSteppingEbos.hpp:769
TimeStepControlType timeStepControl_
time step control object
Definition: AdaptiveTimeSteppingEbos.hpp:757
AdaptiveTimeSteppingEbos(const UnitSystem &unitSystem, const bool terminalOutput=true)
contructor taking parameter object
Definition: AdaptiveTimeSteppingEbos.hpp:265
bool timestepVerbose_
timestep verbosity
Definition: AdaptiveTimeSteppingEbos.hpp:766
bool ignoreConvergenceFailure_
continue instead of stop when minimum time step is reached
Definition: AdaptiveTimeSteppingEbos.hpp:763
double restartFactor_
factor to multiply time step with when solver fails to converge
Definition: AdaptiveTimeSteppingEbos.hpp:758
SimulatorReport step(const SimulatorTimer &simulatorTimer, Solver &solver, const bool isEvent, const std::vector< int > *fipnum=nullptr)
step method that acts like the solver::step method in a sub cycle of time steps
Definition: AdaptiveTimeSteppingEbos.hpp:369
AdaptiveTimeSteppingEbos(const Tuning &tuning, const UnitSystem &unitSystem, const bool terminalOutput=true)
contructor taking parameter object
Definition: AdaptiveTimeSteppingEbos.hpp:292
bool solverVerbose_
solver verbosity
Definition: AdaptiveTimeSteppingEbos.hpp:765
bool useNewtonIteration_
use newton iteration count for adaptive time step control
Definition: AdaptiveTimeSteppingEbos.hpp:770
double suggestedNextTimestep_
suggested size of next timestep
Definition: AdaptiveTimeSteppingEbos.hpp:767
double maxGrowth_
factor that limits the maximum growth of a time step
Definition: AdaptiveTimeSteppingEbos.hpp:760
RelativeChangeInterface.
Definition: TimeStepControlInterface.hpp:32
Definition: SimulatorTimer.hpp:37
double currentStepLength() const override
Current step length.
Definition: SimulatorTimer.cpp:94
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:26
Definition: AdaptiveTimeSteppingEbos.hpp:68
Definition: AdaptiveTimeSteppingEbos.hpp:64
Definition: AdaptiveTimeSteppingEbos.hpp:116
Definition: AdaptiveTimeSteppingEbos.hpp:112
Definition: AdaptiveTimeSteppingEbos.hpp:48
Definition: AdaptiveTimeSteppingEbos.hpp:32
Definition: AdaptiveTimeSteppingEbos.hpp:36
Definition: AdaptiveTimeSteppingEbos.hpp:52
Definition: AdaptiveTimeSteppingEbos.hpp:40
Definition: AdaptiveTimeSteppingEbos.hpp:44
Definition: AdaptiveTimeSteppingEbos.hpp:28
Definition: AdaptiveTimeSteppingEbos.hpp:56
Definition: AdaptiveTimeSteppingEbos.hpp:24
Definition: AdaptiveTimeSteppingEbos.hpp:72
Definition: AdaptiveTimeSteppingEbos.hpp:100
Definition: AdaptiveTimeSteppingEbos.hpp:92
Definition: AdaptiveTimeSteppingEbos.hpp:108
Definition: AdaptiveTimeSteppingEbos.hpp:104
Definition: AdaptiveTimeSteppingEbos.hpp:96
Definition: AdaptiveTimeSteppingEbos.hpp:84
Definition: AdaptiveTimeSteppingEbos.hpp:88
Definition: AdaptiveTimeSteppingEbos.hpp:80
Definition: AdaptiveTimeSteppingEbos.hpp:76
Definition: AdaptiveTimeSteppingEbos.hpp:60
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:31
void reportStep(std::ostringstream &os) const
Print a report suitable for a single simulation step.
Definition: SimulatorReport.cpp:80
Definition: SimulatorReport.hpp:66