3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_TIME_LOOP_HH
25#define DUMUX_TIME_LOOP_HH
26
27#include <algorithm>
28#include <queue>
29#include <iomanip>
30
31#include <dune/common/float_cmp.hh>
32#include <dune/common/timer.hh>
33
34#include <dune/common/version.hh>
35#include <dune/common/parallel/communication.hh>
36
37#include <dune/common/parallel/mpihelper.hh>
38#include <dune/common/exceptions.hh>
39
41#include <dumux/io/format.hh>
42
43namespace Dumux {
44
67template<class Scalar>
69{
70public:
72 virtual ~TimeLoopBase() {};
73
79 virtual Scalar time() const = 0;
80
84 virtual Scalar timeStepSize() const = 0;
85
89 virtual Scalar maxTimeStepSize() const = 0;
90
94 virtual void advanceTimeStep() = 0;
95
100 virtual void setTimeStepSize(Scalar dt) = 0;
101
105 virtual bool finished() const = 0;
106};
107
112template <class Scalar>
113class TimeLoop : public TimeLoopBase<Scalar>
114{
115public:
116 TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
117 : timer_(false)
118 {
119 reset(startTime, dt, tEnd, verbose);
120 }
121
130 void start()
131 {
132 timer_.start();
133 }
134
139 double stop()
140 {
141 return timer_.stop();
142 }
143
148 {
149 timer_.reset();
150 }
151
155 void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
156 {
157 verbose_ =
158 verbose &&
159 Dune::MPIHelper::getCollectiveCommunication().rank() == 0;
160
161 time_ = startTime;
162 endTime_ = tEnd;
163
165 userSetMaxTimeStepSize_ = std::numeric_limits<Scalar>::max();
166 timeStepIdx_ = 0;
167 finished_ = false;
170
171 // ensure that dt is not greater than tEnd-startTime
172 setTimeStepSize(dt);
173
174 timer_.stop();
175 timer_.reset();
176 }
177
181 void advanceTimeStep() override
182 {
183 timeStepIdx_++;
186
187 // compute how long the last time step took
188 const auto cpuTime = wallClockTime();
190 timeAfterLastTimeStep_ = cpuTime;
191
192 // ensure that using current dt we don't exceed tEnd in next time step
194 }
195
202 void setTime(Scalar t)
203 { time_ = t; }
204
211 void setTime(Scalar t, int stepIdx)
212 { time_ = t; timeStepIdx_ = stepIdx; }
213
219 Scalar time() const final
220 { return time_; }
221
225 Scalar endTime() const
226 { return endTime_; }
227
233 void setEndTime(Scalar t)
234 {
235 endTime_ = t;
236 if (verbose_)
237 std::cout << Fmt::format("Set new end time to t = {:.5g} seconds.\n", t);
238 }
239
243 double wallClockTime() const
244 { return timer_.elapsed(); }
245
256 void setTimeStepSize(Scalar dt) final
257 {
258 using std::min;
259 timeStepSize_ = min(dt, maxTimeStepSize());
260 if (!finished() && Dune::FloatCmp::le(timeStepSize_, 0.0, 1e-14*endTime_))
261 std::cerr << Fmt::format("You have set a very small timestep size (dt = {:.5g}).", timeStepSize_)
262 << " This might lead to numerical problems!\n";
263 }
264
271 void setMaxTimeStepSize(Scalar maxDt)
272 {
275 }
276
282 Scalar timeStepSize() const final
283 { return timeStepSize_; }
284
289 int timeStepIndex() const
290 { return timeStepIdx_; }
291
295 Scalar previousTimeStepSize() const
296 { return previousTimeStepSize_; }
297
305 void setFinished(bool finished = true)
306 { finished_ = finished; }
307
314 bool finished() const override
315 {
316 return finished_ || endTime_-time_ < 1e-10*time_;
317 }
318
323 bool willBeFinished() const
324 {
326 }
327
333 Scalar maxTimeStepSize() const override
334 {
335 if (finished())
336 return 0.0;
337
338 using std::min; using std::max;
339 return min(userSetMaxTimeStepSize_, max<Scalar>(0.0, endTime_ - time_));
340 }
341
346 void reportTimeStep() const
347 {
348 if (verbose_)
349 {
350 const auto cpuTime = wallClockTime();
351 using std::round;
352 const auto percent = round( time_ / endTime_ * 100 );
353 std::cout << Fmt::format("[{:3.0f}%] ", percent)
354 << Fmt::format("Time step {} done in {:.2g} seconds. ", timeStepIdx_, timeStepWallClockTime_)
355 << Fmt::format("Wall clock time: {:.2g}, time: {:.5g}, time step size: {:.5g}\n", cpuTime, time_, previousTimeStepSize_);
356 }
357 }
358
362 template< class Communicator = Dune::CollectiveCommunication<typename Dune::MPIHelper::MPICommunicator> >
363 void finalize(const Communicator& comm = Dune::MPIHelper::getCollectiveCommunication())
364 {
365 auto cpuTime = timer_.stop();
366
367 if (verbose_)
368 std::cout << Fmt::format("Simulation took {:.2g} seconds on {} processes.\n", cpuTime, comm.size());
369
370 if (comm.size() > 1)
371 cpuTime = comm.sum(cpuTime);
372
373 if (verbose_)
374 std::cout << Fmt::format("The cumulative CPU time was {:.2g} seconds.\n", cpuTime);
375 }
376
378 bool verbose() const
379 { return verbose_; }
380
382 void setVerbose(bool verbose = true)
383 { verbose_ = verbose; }
384
385 /*
386 * @}
387 */
388
389protected:
390 Dune::Timer timer_;
391 Scalar time_;
392 Scalar endTime_;
393
401};
402
407template <class Scalar>
408class CheckPointTimeLoop : public TimeLoop<Scalar>
409{
410public:
411 CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
412 : TimeLoop<Scalar>(startTime, dt, tEnd, verbose)
413 {
414 periodicCheckPoints_ = false;
415 deltaPeriodicCheckPoint_ = 0.0;
416 lastPeriodicCheckPoint_ = startTime;
417 isCheckPoint_ = false;
418 }
419
423 void advanceTimeStep() override
424 {
425 const auto dt = this->timeStepSize();
426 const auto newTime = this->time()+dt;
427
429 // if we reached a periodic check point
430 if (periodicCheckPoints_ && Dune::FloatCmp::eq(newTime - lastPeriodicCheckPoint_, deltaPeriodicCheckPoint_, 1e-7))
431 {
432 lastPeriodicCheckPoint_ += deltaPeriodicCheckPoint_;
433 isCheckPoint_ = true;
434 }
435
436 // or a manually set check point
437 else if (!checkPoints_.empty() && Dune::FloatCmp::eq(newTime - checkPoints_.front(), 0.0, 1e-7))
438 {
439 checkPoints_.pop();
440 isCheckPoint_ = true;
441 }
442
443 // if not reset the check point flag
444 else
445 {
446 isCheckPoint_ = false;
447 }
448
449 const auto previousTimeStepSize = this->previousTimeStepSize();
450
451 // advance the time step like in the parent class
453
454 // if this is a check point we might have reduced the time step to reach this check point
455 // reset the time step size to the time step size before this time step
456 if (!this->willBeFinished())
457 {
458 using std::max;
459 if (isCheckPoint_)
461
462 // if there is a check point soon check if the time step after the next time step would be smaller
463 // than 20% of the next time step, if yes increase the suggested next time step to exactly reach the check point
464 // (in the limits of the maximum time step size)
465 auto nextDt = this->timeStepSize();
466 const auto threshold = 0.2*nextDt;
467 const auto nextTime = this->time() + nextDt;
468
469 if (periodicCheckPoints_ && Dune::FloatCmp::le(lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - nextTime, threshold))
470 nextDt = lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time();
471
472 if (!checkPoints_.empty() && Dune::FloatCmp::le(checkPoints_.front() - nextTime, threshold))
473 nextDt = checkPoints_.front() - this->time();
474
475 assert(nextDt > 0.0);
476 this->setTimeStepSize(nextDt);
477 }
478 }
479
485 Scalar maxTimeStepSize() const override
486 {
487 using std::min;
488 const auto maxCheckPointDt = computeStepSizeRespectingCheckPoints_();
489 const auto maxDtParent = TimeLoop<Scalar>::maxTimeStepSize();
490 return min(maxDtParent, maxCheckPointDt);
491 }
492
503 void setPeriodicCheckPoint(Scalar interval, Scalar offset = 0.0)
504 {
505 using std::signbit;
506 if (signbit(interval))
507 DUNE_THROW(Dune::InvalidStateException, "Interval has to be positive!");
508
509 periodicCheckPoints_ = true;
510 deltaPeriodicCheckPoint_ = interval;
511 lastPeriodicCheckPoint_ = offset;
512 while (lastPeriodicCheckPoint_ + interval - this->time() < 1e-14*interval)
513 lastPeriodicCheckPoint_ += interval;
514
515 if (this->verbose())
516 std::cout << Fmt::format("Enabled periodic check points every {:.5g} seconds ", interval)
517 << Fmt::format("with the next check point at {:.5g} seconds.\n", lastPeriodicCheckPoint_ + interval);
518
519 // check if the current time point is a check point
520 if (Dune::FloatCmp::eq(this->time()-lastPeriodicCheckPoint_, 0.0, 1e-7))
521 isCheckPoint_ = true;
522
523 // make sure we respect this check point on the next time step
524 this->setTimeStepSize(this->timeStepSize());
525 }
526
529 { periodicCheckPoints_ = false; }
530
533 {
534 periodicCheckPoints_ = false;
535 while (!checkPoints_.empty())
536 checkPoints_.pop();
537 }
538
543 bool isCheckPoint() const
544 { return isCheckPoint_; }
545
552 void setCheckPoint(Scalar t)
553 {
554 // set the check point
555 setCheckPoint_(t);
556
557 // make sure we respect this check point on the next time step
558 this->setTimeStepSize(this->timeStepSize());
559 }
560
567 void setCheckPoint(const std::vector<Scalar>& checkPoints)
568 { setCheckPoint(checkPoints.begin(), checkPoints.end()); }
569
577 template<class ForwardIterator>
578 void setCheckPoint(ForwardIterator first, ForwardIterator last)
579 {
580 // set the check points
581 for (; first != last; ++first)
582 setCheckPoint_(*first);
583
584 // make sure we respect this check point on the next time step
585 this->setTimeStepSize(this->timeStepSize());
586 }
587
588private:
590 void setCheckPoint_(Scalar t)
591 {
592 if (Dune::FloatCmp::le(t - this->time(), 0.0, this->timeStepSize()*1e-7))
593 {
594 if (this->verbose())
595 std::cerr << Fmt::format("Couldn't insert checkpoint at t = {:.5g} ", t)
596 << Fmt::format("because that's in the past! (current simulation time is {:.5g})\n", this->time());
597 return;
598 }
599
600 if (!checkPoints_.empty())
601 {
602 if (t < checkPoints_.back())
603 {
604 if (this->verbose())
605 std::cerr << "Couldn't insert checkpoint as it is earlier than the last check point in the queue.\n"
606 << "Checkpoints can only be inserted in ascending order." << std::endl;
607 return;
608 }
609 }
610
611 checkPoints_.push(t);
612 if (this->verbose())
613 std::cout << Fmt::format("Set check point at t = {:.5g} seconds.\n", t);
614 }
615
619 Scalar computeStepSizeRespectingCheckPoints_() const
620 {
621 using std::min;
622 auto maxDt = std::numeric_limits<Scalar>::max();
623
624 if (periodicCheckPoints_)
625 maxDt = min(maxDt, lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time());
626
627 if (!checkPoints_.empty())
628 maxDt = min(maxDt, checkPoints_.front() - this->time());
629
630 return maxDt;
631 }
632
633 bool periodicCheckPoints_;
634 Scalar deltaPeriodicCheckPoint_;
635 Scalar lastPeriodicCheckPoint_;
636 std::queue<Scalar> checkPoints_;
637 bool isCheckPoint_;
638};
639
640} // end namespace Dumux
641
642#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Formatting based on the fmt-library which implements std::format of C++20.
Definition: adapt.hh:29
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
virtual ~TimeLoopBase()
Abstract base class needs virtual constructor.
Definition: common/timeloop.hh:72
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:114
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: common/timeloop.hh:333
void setMaxTimeStepSize(Scalar maxDt)
Set the maximum time step size to a given value.
Definition: common/timeloop.hh:271
void resetTimer()
Reset the timer.
Definition: common/timeloop.hh:147
TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: common/timeloop.hh:116
void setVerbose(bool verbose=true)
Sets time loop verbosity.
Definition: common/timeloop.hh:382
void setTime(Scalar t, int stepIdx)
Set the current simulated time and the time step index.
Definition: common/timeloop.hh:211
Scalar time_
Definition: common/timeloop.hh:391
int timeStepIndex() const
Returns number of time steps which have been executed since the beginning of the simulation.
Definition: common/timeloop.hh:289
double stop()
Tells the time loop to stop tracking the time.
Definition: common/timeloop.hh:139
Scalar endTime_
Definition: common/timeloop.hh:392
Scalar previousTimeStepSize_
Definition: common/timeloop.hh:395
void advanceTimeStep() override
Advance time step.
Definition: common/timeloop.hh:181
int timeStepIdx_
Definition: common/timeloop.hh:398
void start()
Tells the time loop to start tracking the time.
Definition: common/timeloop.hh:130
Scalar userSetMaxTimeStepSize_
Definition: common/timeloop.hh:396
void setFinished(bool finished=true)
Specify whether the simulation is finished.
Definition: common/timeloop.hh:305
bool verbose() const
If the time loop has verbose output.
Definition: common/timeloop.hh:378
bool verbose_
Definition: common/timeloop.hh:400
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:323
Scalar previousTimeStepSize() const
The previous time step size.
Definition: common/timeloop.hh:295
void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Reset the time loop.
Definition: common/timeloop.hh:155
Scalar endTime() const
Returns the number of (simulated) seconds which the simulation runs.
Definition: common/timeloop.hh:225
void setTimeStepSize(Scalar dt) final
Set the current time step size to a given value.
Definition: common/timeloop.hh:256
bool finished() const override
Returns true if the simulation is finished.
Definition: common/timeloop.hh:314
Scalar timeStepSize_
Definition: common/timeloop.hh:394
void finalize(const Communicator &comm=Dune::MPIHelper::getCollectiveCommunication())
Print final status and stops tracking the time.
Definition: common/timeloop.hh:363
void reportTimeStep() const
State info on cpu time.
Definition: common/timeloop.hh:346
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:282
bool finished_
Definition: common/timeloop.hh:399
Dune::Timer timer_
Definition: common/timeloop.hh:390
void setTime(Scalar t)
Set the current simulated time, don't change the current time step index.
Definition: common/timeloop.hh:202
void setEndTime(Scalar t)
Set the time of simulated seconds at which the simulation runs.
Definition: common/timeloop.hh:233
double wallClockTime() const
Returns the current wall clock time (cpu time) spend in this time loop.
Definition: common/timeloop.hh:243
Scalar timeAfterLastTimeStep_
Definition: common/timeloop.hh:397
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:219
Scalar timeStepWallClockTime_
Definition: common/timeloop.hh:397
A time loop with a check point mechanism.
Definition: common/timeloop.hh:409
bool isCheckPoint() const
Whether now is a time checkpoint.
Definition: common/timeloop.hh:543
void setPeriodicCheckPoint(Scalar interval, Scalar offset=0.0)
Set a periodic check point.
Definition: common/timeloop.hh:503
void setCheckPoint(Scalar t)
add a checkpoint to the queue
Definition: common/timeloop.hh:552
void advanceTimeStep() override
Advance time step.
Definition: common/timeloop.hh:423
void removeAllCheckPoints()
remove all check points
Definition: common/timeloop.hh:532
void setCheckPoint(const std::vector< Scalar > &checkPoints)
add checkpoints to the queue from a vector of time points
Definition: common/timeloop.hh:567
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: common/timeloop.hh:485
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:578
CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: common/timeloop.hh:411
void disablePeriodicCheckPoints()
disable periodic check points
Definition: common/timeloop.hh:528