version 3.7
common/timeloop.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_TIME_LOOP_HH
13#define DUMUX_TIME_LOOP_HH
14
15#include <algorithm>
16#include <queue>
17#include <iomanip>
18
19#include <dune/common/float_cmp.hh>
20#include <dune/common/timer.hh>
21
22#include <dune/common/parallel/communication.hh>
23
24#include <dune/common/parallel/mpihelper.hh>
25#include <dune/common/exceptions.hh>
26
28#include <dumux/io/format.hh>
29
30namespace Dumux {
31
54template<class Scalar>
56{
57public:
59 virtual ~TimeLoopBase() {};
60
66 virtual Scalar time() const = 0;
67
71 virtual Scalar timeStepSize() const = 0;
72
76 virtual Scalar maxTimeStepSize() const = 0;
77
81 virtual void advanceTimeStep() = 0;
82
87 virtual void setTimeStepSize(Scalar dt) = 0;
88
92 virtual bool finished() const = 0;
93};
94
99template <class Scalar>
100class TimeLoop : public TimeLoopBase<Scalar>
101{
102public:
103 TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
104 : timer_(false)
105 {
106 reset(startTime, dt, tEnd, verbose);
107 }
108
117 void start()
118 {
119 timer_.start();
120 }
121
126 double stop()
127 {
128 return timer_.stop();
129 }
130
135 {
136 timer_.reset();
137 }
138
142 void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
143 {
144 verbose_ =
145 verbose &&
146 Dune::MPIHelper::getCommunication().rank() == 0;
147
148 time_ = startTime;
149 endTime_ = tEnd;
150
152 userSetMaxTimeStepSize_ = std::numeric_limits<Scalar>::max();
153 timeStepIdx_ = 0;
154 finished_ = false;
157
158 // ensure that dt is not greater than tEnd-startTime
159 setTimeStepSize(dt);
160
161 timer_.stop();
162 timer_.reset();
163 }
164
168 void advanceTimeStep() override
169 {
170 timeStepIdx_++;
173
174 // compute how long the last time step took
175 const auto cpuTime = wallClockTime();
177 timeAfterLastTimeStep_ = cpuTime;
178
179 // ensure that using current dt we don't exceed tEnd in next time step
181 }
182
189 void setTime(Scalar t)
190 { time_ = t; }
191
198 void setTime(Scalar t, int stepIdx)
199 { time_ = t; timeStepIdx_ = stepIdx; }
200
206 Scalar time() const final
207 { return time_; }
208
212 Scalar endTime() const
213 { return endTime_; }
214
220 void setEndTime(Scalar t)
221 {
222 endTime_ = t;
223 if (verbose_)
224 std::cout << Fmt::format("Set new end time to t = {:.5g} seconds.\n", t);
225 }
226
230 double wallClockTime() const
231 { return timer_.elapsed(); }
232
243 void setTimeStepSize(Scalar dt) final
244 {
245 using std::min;
246 timeStepSize_ = min(dt, maxTimeStepSize());
247 if (!finished() && Dune::FloatCmp::le(timeStepSize_, 0.0, 1e-14*endTime_))
248 std::cerr << Fmt::format("You have set a very small timestep size (dt = {:.5g}).", timeStepSize_)
249 << " This might lead to numerical problems!\n";
250 }
251
258 void setMaxTimeStepSize(Scalar maxDt)
259 {
262 }
263
269 Scalar timeStepSize() const final
270 { return timeStepSize_; }
271
276 int timeStepIndex() const
277 { return timeStepIdx_; }
278
282 Scalar previousTimeStepSize() const
283 { return previousTimeStepSize_; }
284
292 void setFinished(bool finished = true)
293 { finished_ = finished; }
294
301 bool finished() const override
302 {
303 return finished_ || endTime_-time_ < 1e-10*time_;
304 }
305
310 bool willBeFinished() const
311 {
313 }
314
320 Scalar maxTimeStepSize() const override
321 {
322 if (finished())
323 return 0.0;
324
325 using std::min; using std::max;
326 return min(userSetMaxTimeStepSize_, max<Scalar>(0.0, endTime_ - time_));
327 }
328
333 void reportTimeStep() const
334 {
335 if (verbose_)
336 {
337 const auto cpuTime = wallClockTime();
338 using std::round;
339 const auto percent = round( time_ / endTime_ * 100 );
340 std::cout << Fmt::format("[{:3.0f}%] ", percent)
341 << Fmt::format("Time step {} done in {:.2g} seconds. ", timeStepIdx_, timeStepWallClockTime_)
342 << Fmt::format("Wall clock time: {:.5g}, time: {:.5g}, time step size: {:.5g}\n", cpuTime, time_, previousTimeStepSize_);
343 }
344 }
345
349 template< class Communicator = Dune::Communication<typename Dune::MPIHelper::MPICommunicator> >
350 void finalize(const Communicator& comm = Dune::MPIHelper::getCommunication())
351 {
352 auto cpuTime = timer_.stop();
353
354 if (verbose_)
355 std::cout << Fmt::format("Simulation took {:.5g} seconds on {} processes.\n", cpuTime, comm.size());
356
357 if (comm.size() > 1)
358 cpuTime = comm.sum(cpuTime);
359
360 if (verbose_)
361 std::cout << Fmt::format("The cumulative CPU time was {:.5g} seconds.\n", cpuTime);
362 }
363
365 bool verbose() const
366 { return verbose_; }
367
369 void setVerbose(bool verbose = true)
370 { verbose_ = verbose; }
371
372 /*
373 * @}
374 */
375
376protected:
377 Dune::Timer timer_;
378 Scalar time_;
379 Scalar endTime_;
380
388};
389
394template <class Scalar>
395class CheckPointTimeLoop : public TimeLoop<Scalar>
396{
397public:
398 CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
399 : TimeLoop<Scalar>(startTime, dt, tEnd, verbose)
400 {
401 periodicCheckPoints_ = false;
402 deltaPeriodicCheckPoint_ = 0.0;
403 lastPeriodicCheckPoint_ = startTime;
404 isCheckPoint_ = false;
405 }
406
410 void advanceTimeStep() override
411 {
412 const auto dt = this->timeStepSize();
413 const auto newTime = this->time()+dt;
414
416 // if we reached a periodic check point
417 if (periodicCheckPoints_ && Dune::FloatCmp::eq(newTime - lastPeriodicCheckPoint_, deltaPeriodicCheckPoint_, 1e-7))
418 {
419 lastPeriodicCheckPoint_ += deltaPeriodicCheckPoint_;
420 isCheckPoint_ = true;
421 }
422
423 // or a manually set check point
424 else if (!checkPoints_.empty() && Dune::FloatCmp::eq(newTime - checkPoints_.front(), 0.0, 1e-7))
425 {
426 checkPoints_.pop();
427 isCheckPoint_ = true;
428 }
429
430 // if not reset the check point flag
431 else
432 {
433 isCheckPoint_ = false;
434 }
435
436 const auto previousTimeStepSize = this->previousTimeStepSize();
437
438 // advance the time step like in the parent class
440
441 // if this is a check point we might have reduced the time step to reach this check point
442 // reset the time step size to the time step size before this time step
443 if (!this->willBeFinished())
444 {
445 using std::max;
446 if (isCheckPoint_)
448
449 // if there is a check point soon check if the time step after the next time step would be smaller
450 // than 20% of the next time step, if yes increase the suggested next time step to exactly reach the check point
451 // (in the limits of the maximum time step size)
452 auto nextDt = this->timeStepSize();
453 const auto threshold = 0.2*nextDt;
454 const auto nextTime = this->time() + nextDt;
455
456 if (periodicCheckPoints_ && Dune::FloatCmp::le(lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - nextTime, threshold))
457 nextDt = lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time();
458
459 if (!checkPoints_.empty() && Dune::FloatCmp::le(checkPoints_.front() - nextTime, threshold))
460 nextDt = checkPoints_.front() - this->time();
461
462 assert(nextDt > 0.0);
463 this->setTimeStepSize(nextDt);
464 }
465 }
466
472 Scalar maxTimeStepSize() const override
473 {
474 using std::min;
475 const auto maxCheckPointDt = computeStepSizeRespectingCheckPoints_();
476 const auto maxDtParent = TimeLoop<Scalar>::maxTimeStepSize();
477 return min(maxDtParent, maxCheckPointDt);
478 }
479
490 void setPeriodicCheckPoint(Scalar interval, Scalar offset = 0.0)
491 {
492 using std::signbit;
493 if (signbit(interval))
494 DUNE_THROW(Dune::InvalidStateException, "Interval has to be positive!");
495
496 periodicCheckPoints_ = true;
497 deltaPeriodicCheckPoint_ = interval;
498 lastPeriodicCheckPoint_ = offset;
499 while (lastPeriodicCheckPoint_ + interval - this->time() < 1e-14*interval)
500 lastPeriodicCheckPoint_ += interval;
501
502 if (this->verbose())
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);
505
506 // check if the current time point is a check point
507 if (Dune::FloatCmp::eq(this->time()-lastPeriodicCheckPoint_, 0.0, 1e-7))
508 isCheckPoint_ = true;
509
510 // make sure we respect this check point on the next time step
511 this->setTimeStepSize(this->timeStepSize());
512 }
513
516 { periodicCheckPoints_ = false; }
517
520 {
521 periodicCheckPoints_ = false;
522 while (!checkPoints_.empty())
523 checkPoints_.pop();
524 }
525
530 bool isCheckPoint() const
531 { return isCheckPoint_; }
532
539 void setCheckPoint(Scalar t)
540 {
541 // set the check point
542 setCheckPoint_(t);
543
544 // make sure we respect this check point on the next time step
545 this->setTimeStepSize(this->timeStepSize());
546 }
547
554 void setCheckPoint(const std::vector<Scalar>& checkPoints)
555 { setCheckPoint(checkPoints.begin(), checkPoints.end()); }
556
564 template<class ForwardIterator>
565 void setCheckPoint(ForwardIterator first, ForwardIterator last)
566 {
567 // set the check points
568 for (; first != last; ++first)
569 setCheckPoint_(*first);
570
571 // make sure we respect this check point on the next time step
572 this->setTimeStepSize(this->timeStepSize());
573 }
574
575private:
577 void setCheckPoint_(Scalar t)
578 {
579 if (Dune::FloatCmp::le(t - this->time(), 0.0, this->timeStepSize()*1e-7))
580 {
581 if (this->verbose())
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());
584 return;
585 }
586
587 if (!checkPoints_.empty())
588 {
589 if (t <= checkPoints_.back())
590 {
591 if (this->verbose())
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;
595 return;
596 }
597 }
598
599 checkPoints_.push(t);
600 if (this->verbose())
601 std::cout << Fmt::format("Set check point at t = {:.5g} seconds.\n", t);
602 }
603
607 Scalar computeStepSizeRespectingCheckPoints_() const
608 {
609 using std::min;
610 auto maxDt = std::numeric_limits<Scalar>::max();
611
612 if (periodicCheckPoints_)
613 maxDt = min(maxDt, lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time());
614
615 if (!checkPoints_.empty())
616 maxDt = min(maxDt, checkPoints_.front() - this->time());
617
618 return maxDt;
619 }
620
621 bool periodicCheckPoints_;
622 Scalar deltaPeriodicCheckPoint_;
623 Scalar lastPeriodicCheckPoint_;
624 std::queue<Scalar> checkPoints_;
625 bool isCheckPoint_;
626};
627
628} // end namespace Dumux
629
630#endif
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
Formatting based on the fmt-library which implements std::format of C++20.
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.