12#ifndef DUMUX_TIME_LOOP_HH
13#define DUMUX_TIME_LOOP_HH
19#include <dune/common/float_cmp.hh>
20#include <dune/common/timer.hh>
22#include <dune/common/parallel/communication.hh>
24#include <dune/common/parallel/mpihelper.hh>
25#include <dune/common/exceptions.hh>
66 virtual Scalar
time()
const = 0;
99template <
class Scalar>
142 void reset(Scalar startTime, Scalar dt, Scalar tEnd,
bool verbose =
true)
146 Dune::MPIHelper::getCommunication().rank() == 0;
224 std::cout << Fmt::format(
"Set new end time to t = {:.5g} seconds.\n", t);
231 {
return timer_.elapsed(); }
248 std::cerr << Fmt::format(
"You have set a very small timestep size (dt = {:.5g}).",
timeStepSize_)
249 <<
" This might lead to numerical problems!\n";
325 using std::min;
using std::max;
340 std::cout << Fmt::format(
"[{:3.0f}%] ", percent)
349 template<
class Communicator = Dune::Communication<
typename Dune::MPIHelper::MPICommunicator> >
350 void finalize(
const Communicator& comm = Dune::MPIHelper::getCommunication())
352 auto cpuTime =
timer_.stop();
355 std::cout << Fmt::format(
"Simulation took {:.5g} seconds on {} processes.\n", cpuTime, comm.size());
358 cpuTime = comm.sum(cpuTime);
361 std::cout << Fmt::format(
"The cumulative CPU time was {:.5g} seconds.\n", cpuTime);
394template <
class Scalar>
401 periodicCheckPoints_ =
false;
402 deltaPeriodicCheckPoint_ = 0.0;
403 lastPeriodicCheckPoint_ = startTime;
404 isCheckPoint_ =
false;
413 const auto newTime = this->
time()+dt;
417 if (periodicCheckPoints_ && Dune::FloatCmp::eq(newTime - lastPeriodicCheckPoint_, deltaPeriodicCheckPoint_, 1e-7))
419 lastPeriodicCheckPoint_ += deltaPeriodicCheckPoint_;
420 isCheckPoint_ =
true;
424 else if (!checkPoints_.empty() && Dune::FloatCmp::eq(newTime - checkPoints_.front(), 0.0, 1e-7))
427 isCheckPoint_ =
true;
433 isCheckPoint_ =
false;
453 const auto threshold = 0.2*nextDt;
454 const auto nextTime = this->
time() + nextDt;
456 if (periodicCheckPoints_ && Dune::FloatCmp::le(lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - nextTime, threshold))
457 nextDt = lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->
time();
459 if (!checkPoints_.empty() && Dune::FloatCmp::le(checkPoints_.front() - nextTime, threshold))
460 nextDt = checkPoints_.front() - this->
time();
462 assert(nextDt > 0.0);
475 const auto maxCheckPointDt = computeStepSizeRespectingCheckPoints_();
477 return min(maxDtParent, maxCheckPointDt);
493 if (signbit(interval))
494 DUNE_THROW(Dune::InvalidStateException,
"Interval has to be positive!");
496 periodicCheckPoints_ =
true;
497 deltaPeriodicCheckPoint_ = interval;
498 lastPeriodicCheckPoint_ = offset;
499 while (lastPeriodicCheckPoint_ + interval - this->
time() < 1e-14*interval)
500 lastPeriodicCheckPoint_ += interval;
503 std::cout << Fmt::format(
"Enabled periodic check points every {:.5g} seconds ", interval)
504 << Fmt::format(
"with the next check point at {:.5g} seconds.\n", lastPeriodicCheckPoint_ + interval);
507 if (Dune::FloatCmp::eq(this->
time()-lastPeriodicCheckPoint_, 0.0, 1e-7))
508 isCheckPoint_ =
true;
516 { periodicCheckPoints_ =
false; }
521 periodicCheckPoints_ =
false;
522 while (!checkPoints_.empty())
531 {
return isCheckPoint_; }
564 template<
class ForwardIterator>
568 for (; first != last; ++first)
569 setCheckPoint_(*first);
577 void setCheckPoint_(Scalar t)
582 std::cerr << Fmt::format(
"Couldn't insert checkpoint at t = {:.5g} ", t)
583 << Fmt::format(
"because that's in the past! (current simulation time is {:.5g})\n", this->
time());
587 if (!checkPoints_.empty())
589 if (t <= checkPoints_.back())
592 std::cerr << Fmt::format(
"Couldn't insert checkpoint at t = {:.5g} ", t)
593 << Fmt::format(
"because it's earlier than or equal to the last check point (t = {:.5g}) in the queue.\n", checkPoints_.back())
594 <<
"Checkpoints can only be inserted in ascending order." << std::endl;
599 checkPoints_.push(t);
601 std::cout << Fmt::format(
"Set check point at t = {:.5g} seconds.\n", t);
607 Scalar computeStepSizeRespectingCheckPoints_()
const
610 auto maxDt = std::numeric_limits<Scalar>::max();
612 if (periodicCheckPoints_)
613 maxDt = min(maxDt, lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->
time());
615 if (!checkPoints_.empty())
616 maxDt = min(maxDt, checkPoints_.front() - this->time());
621 bool periodicCheckPoints_;
622 Scalar deltaPeriodicCheckPoint_;
623 Scalar lastPeriodicCheckPoint_;
624 std::queue<Scalar> checkPoints_;
A time loop with a check point mechanism.
Definition: common/timeloop.hh:396
bool isCheckPoint() const
Whether now is a time checkpoint.
Definition: common/timeloop.hh:530
void setPeriodicCheckPoint(Scalar interval, Scalar offset=0.0)
Set a periodic check point.
Definition: common/timeloop.hh:490
void setCheckPoint(Scalar t)
add a checkpoint to the queue
Definition: common/timeloop.hh:539
void advanceTimeStep() override
Advance time step.
Definition: common/timeloop.hh:410
void removeAllCheckPoints()
remove all check points
Definition: common/timeloop.hh:519
void setCheckPoint(const std::vector< Scalar > &checkPoints)
add checkpoints to the queue from a vector of time points
Definition: common/timeloop.hh:554
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: common/timeloop.hh:472
void setCheckPoint(ForwardIterator first, ForwardIterator last)
add checkpoints to the queue from a container from the first iterator to the last iterator
Definition: common/timeloop.hh:565
CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: common/timeloop.hh:398
void disablePeriodicCheckPoints()
disable periodic check points
Definition: common/timeloop.hh:515
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:56
virtual ~TimeLoopBase()
Abstract base class needs virtual constructor.
Definition: common/timeloop.hh:59
virtual void setTimeStepSize(Scalar dt)=0
Set the current time step size to a given value.
virtual Scalar time() const =0
Return the time before the time integration. To get the time after the time integration you have to ...
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
virtual bool finished() const =0
Returns true if the simulation is finished.
virtual void advanceTimeStep()=0
Advance to the next time step.
virtual Scalar maxTimeStepSize() const =0
Get the maximum possible time step size .
The default time loop for instationary simulations.
Definition: common/timeloop.hh:101
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: common/timeloop.hh:320
void setMaxTimeStepSize(Scalar maxDt)
Set the maximum time step size to a given value.
Definition: common/timeloop.hh:258
void resetTimer()
Reset the timer.
Definition: common/timeloop.hh:134
TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: common/timeloop.hh:103
void setVerbose(bool verbose=true)
Sets time loop verbosity.
Definition: common/timeloop.hh:369
void setTime(Scalar t, int stepIdx)
Set the current simulated time and the time step index.
Definition: common/timeloop.hh:198
Scalar time_
Definition: common/timeloop.hh:378
int timeStepIndex() const
Returns number of time steps which have been executed since the beginning of the simulation.
Definition: common/timeloop.hh:276
double stop()
Tells the time loop to stop tracking the time.
Definition: common/timeloop.hh:126
Scalar endTime_
Definition: common/timeloop.hh:379
Scalar previousTimeStepSize_
Definition: common/timeloop.hh:382
void advanceTimeStep() override
Advance time step.
Definition: common/timeloop.hh:168
int timeStepIdx_
Definition: common/timeloop.hh:385
void start()
Tells the time loop to start tracking the time.
Definition: common/timeloop.hh:117
Scalar userSetMaxTimeStepSize_
Definition: common/timeloop.hh:383
void setFinished(bool finished=true)
Specify whether the simulation is finished.
Definition: common/timeloop.hh:292
bool verbose() const
If the time loop has verbose output.
Definition: common/timeloop.hh:365
void finalize(const Communicator &comm=Dune::MPIHelper::getCommunication())
Print final status and stops tracking the time.
Definition: common/timeloop.hh:350
bool verbose_
Definition: common/timeloop.hh:387
bool willBeFinished() const
Returns true if the simulation is finished after the time level is incremented by the current time st...
Definition: common/timeloop.hh:310
Scalar previousTimeStepSize() const
The previous time step size.
Definition: common/timeloop.hh:282
void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Reset the time loop.
Definition: common/timeloop.hh:142
Scalar endTime() const
Returns the number of (simulated) seconds which the simulation runs.
Definition: common/timeloop.hh:212
void setTimeStepSize(Scalar dt) final
Set the current time step size to a given value.
Definition: common/timeloop.hh:243
bool finished() const override
Returns true if the simulation is finished.
Definition: common/timeloop.hh:301
Scalar timeStepSize_
Definition: common/timeloop.hh:381
void reportTimeStep() const
State info on cpu time.
Definition: common/timeloop.hh:333
Scalar timeStepSize() const final
Returns the suggested time step length so that we don't miss the beginning of the next episode or cr...
Definition: common/timeloop.hh:269
bool finished_
Definition: common/timeloop.hh:386
Dune::Timer timer_
Definition: common/timeloop.hh:377
void setTime(Scalar t)
Set the current simulated time, don't change the current time step index.
Definition: common/timeloop.hh:189
void setEndTime(Scalar t)
Set the time of simulated seconds at which the simulation runs.
Definition: common/timeloop.hh:220
double wallClockTime() const
Returns the current wall clock time (cpu time) spend in this time loop.
Definition: common/timeloop.hh:230
Scalar timeAfterLastTimeStep_
Definition: common/timeloop.hh:384
Scalar time() const final
Return the time before the time integration. To get the time after the time integration you have to ...
Definition: common/timeloop.hh:206
Scalar timeStepWallClockTime_
Definition: common/timeloop.hh:384
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.