3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
36#include <dune/common/parallel/collectivecommunication.hh>
37#else
38#include <dune/common/parallel/communication.hh>
39#endif
40
41#include <dune/common/parallel/mpihelper.hh>
42#include <dune/common/exceptions.hh>
43
45
46namespace Dumux {
47
70template<class Scalar>
72{
73public:
75 virtual ~TimeLoopBase() {};
76
82 virtual Scalar time() const = 0;
83
87 virtual Scalar timeStepSize() const = 0;
88
92 virtual Scalar maxTimeStepSize() const = 0;
93
97 virtual void advanceTimeStep() = 0;
98
103 virtual void setTimeStepSize(Scalar dt) = 0;
104
108 virtual bool finished() const = 0;
109};
110
115template <class Scalar>
116class TimeLoop : public TimeLoopBase<Scalar>
117{
118public:
119 TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
120 : timer_(false)
121 {
122 reset(startTime, dt, tEnd, verbose);
123 }
124
133 void start()
134 {
135 timer_.start();
136 }
137
142 double stop()
143 {
144 return timer_.stop();
145 }
146
151 {
152 timer_.reset();
153 }
154
158 void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
159 {
160 verbose_ =
161 verbose &&
162 Dune::MPIHelper::getCollectiveCommunication().rank() == 0;
163
164 time_ = startTime;
165 endTime_ = tEnd;
166
168 userSetMaxTimeStepSize_ = std::numeric_limits<Scalar>::max();
169 timeStepIdx_ = 0;
170 finished_ = false;
173
174 // ensure that dt is not greater than tEnd-startTime
175 setTimeStepSize(dt);
176
177 timer_.stop();
178 timer_.reset();
179 }
180
184 void advanceTimeStep() override
185 {
186 timeStepIdx_++;
189
190 // compute how long the last time step took
191 const auto cpuTime = wallClockTime();
193 timeAfterLastTimeStep_ = cpuTime;
194
195 // ensure that using current dt we don't exceed tEnd in next time step
197 }
198
205 void setTime(Scalar t)
206 { time_ = t; }
207
214 void setTime(Scalar t, int stepIdx)
215 { time_ = t; timeStepIdx_ = stepIdx; }
216
222 Scalar time() const final
223 { return time_; }
224
228 Scalar endTime() const
229 { return endTime_; }
230
236 void setEndTime(Scalar t)
237 {
238 endTime_ = t;
239 if (verbose_)
240 std::cout << "Set new end time to t = " << t << " seconds." << std::endl;
241 }
242
246 double wallClockTime() const
247 { return timer_.elapsed(); }
248
259 void setTimeStepSize(Scalar dt) final
260 {
261 using std::min;
262 timeStepSize_ = min(dt, maxTimeStepSize());
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 const auto percent = std::round( time_ / endTime_ * 100 );
352 std::cout << "[" << std::fixed << std::setw( 3 ) << std::setfill( ' ' )
353 << std::setprecision( 0 ) << percent << "%] "
354 << "Time step " << timeStepIdx_ << " done in "
355 << std::setprecision( 6 ) << timeStepWallClockTime_ << " seconds. "
356 << "Wall clock time: " << std::setprecision( 3 ) << cpuTime
357 << ", time: " << std::setprecision( 5 ) << time_
358 << ", time step size: " << std::setprecision( 8 ) << previousTimeStepSize_
359 << std::endl;
360 }
361 }
362
366 template< class Communicator = Dune::CollectiveCommunication<typename Dune::MPIHelper::MPICommunicator> >
367 void finalize(const Communicator& comm = Dune::MPIHelper::getCollectiveCommunication())
368 {
369 auto cpuTime = timer_.stop();
370
371 if (verbose_)
372 {
373 std::cout << "Simulation took " << cpuTime << " seconds on "
374 << comm.size() << " processes.\n";
375 }
376
377 if (comm.size() > 1)
378 cpuTime = comm.sum(cpuTime);
379
380 if (verbose_)
381 {
382 std::cout << "The cumulative CPU time was " << cpuTime << " seconds.\n";
383 }
384 }
385
387 bool verbose() const
388 { return verbose_; }
389
390 /*
391 * @}
392 */
393
394protected:
395 Dune::Timer timer_;
396 Scalar time_;
397 Scalar endTime_;
398
406};
407
412template <class Scalar>
413class CheckPointTimeLoop : public TimeLoop<Scalar>
414{
415public:
416 CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
417 : TimeLoop<Scalar>(startTime, dt, tEnd, verbose)
418 {
419 periodicCheckPoints_ = false;
420 deltaPeriodicCheckPoint_ = 0.0;
421 lastPeriodicCheckPoint_ = startTime;
422 isCheckPoint_ = false;
423 }
424
428 void advanceTimeStep() override
429 {
430 const auto dt = this->timeStepSize();
431 const auto newTime = this->time()+dt;
432
434 // if we reached a periodic check point
435 if (periodicCheckPoints_ && Dune::FloatCmp::eq(newTime - lastPeriodicCheckPoint_, deltaPeriodicCheckPoint_, 1e-7))
436 {
437 lastPeriodicCheckPoint_ += deltaPeriodicCheckPoint_;
438 isCheckPoint_ = true;
439 }
440
441 // or a manually set check point
442 else if (!checkPoints_.empty() && Dune::FloatCmp::eq(newTime - checkPoints_.front(), 0.0, 1e-7))
443 {
444 checkPoints_.pop();
445 isCheckPoint_ = true;
446 }
447
448 // if not reset the check point flag
449 else
450 {
451 isCheckPoint_ = false;
452 }
453
454 const auto previousTimeStepSize = this->previousTimeStepSize();
455
456 // advance the time step like in the parent class
458
459 // if this is a check point we might have reduced the time step to reach this check point
460 // reset the time step size to the time step size before this time step
461 if (!this->willBeFinished())
462 {
463 using std::max;
464 if (isCheckPoint_)
466
467 // if there is a check point soon check if the time step after the next time step would be smaller
468 // than 20% of the next time step, if yes increase the suggested next time step to exactly reach the check point
469 // (in the limits of the maximum time step size)
470 auto nextDt = this->timeStepSize();
471 const auto threshold = 0.2*nextDt;
472 const auto nextTime = this->time() + nextDt;
473
474 if (periodicCheckPoints_ && Dune::FloatCmp::le(lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - nextTime, threshold))
475 nextDt = lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time();
476
477 if (!checkPoints_.empty() && Dune::FloatCmp::le(checkPoints_.front() - nextTime, threshold))
478 nextDt = checkPoints_.front() - this->time();
479
480 assert(nextDt > 0.0);
481 this->setTimeStepSize(nextDt);
482 }
483 }
484
490 Scalar maxTimeStepSize() const override
491 {
492 using std::min;
493 const auto maxCheckPointDt = computeStepSizeRespectingCheckPoints_();
494 const auto maxDtParent = TimeLoop<Scalar>::maxTimeStepSize();
495 return min(maxDtParent, maxCheckPointDt);
496 }
497
508 void setPeriodicCheckPoint(Scalar interval, Scalar offset = 0.0)
509 {
510 using std::signbit;
511 if (signbit(interval))
512 DUNE_THROW(Dune::InvalidStateException, "Interval has to be positive!");
513
514 periodicCheckPoints_ = true;
515 deltaPeriodicCheckPoint_ = interval;
516 lastPeriodicCheckPoint_ = offset;
517 while (lastPeriodicCheckPoint_ + interval - this->time() < 1e-14*interval)
518 lastPeriodicCheckPoint_ += interval;
519
520 if (this->verbose())
521 std::cout << "Enabled periodic check points every " << interval
522 << " seconds with the next check point at " << lastPeriodicCheckPoint_ + interval << " seconds." << std::endl;
523
524 // check if the current time point is a check point
525 if (Dune::FloatCmp::eq(this->time()-lastPeriodicCheckPoint_, 0.0, 1e-7))
526 isCheckPoint_ = true;
527
528 // make sure we respect this check point on the next time step
529 this->setTimeStepSize(this->timeStepSize());
530 }
531
534 { periodicCheckPoints_ = false; }
535
538 {
539 periodicCheckPoints_ = false;
540 while (!checkPoints_.empty())
541 checkPoints_.pop();
542 }
543
548 bool isCheckPoint() const
549 { return isCheckPoint_; }
550
557 void setCheckPoint(Scalar t)
558 {
559 // set the check point
560 setCheckPoint_(t);
561
562 // make sure we respect this check point on the next time step
563 this->setTimeStepSize(this->timeStepSize());
564 }
565
572 void setCheckPoint(const std::vector<Scalar>& checkPoints)
573 { setCheckPoint(checkPoints.begin(), checkPoints.end()); }
574
582 template<class ForwardIterator>
583 void setCheckPoint(ForwardIterator first, ForwardIterator last)
584 {
585 // set the check points
586 for (; first != last; ++first)
587 setCheckPoint_(*first);
588
589 // make sure we respect this check point on the next time step
590 this->setTimeStepSize(this->timeStepSize());
591 }
592
593private:
595 void setCheckPoint_(Scalar t)
596 {
597 if (t < this->time())
598 {
599 if (this->verbose())
600 std::cerr << "Couldn't insert checkpoint at t = " << t
601 << " because that's in the past! (current simulation time is " << this->time() << ")" << std::endl;
602 return;
603 }
604
605 if (!checkPoints_.empty())
606 {
607 if (t < checkPoints_.back())
608 {
609 if (this->verbose())
610 std::cerr << "Couldn't insert checkpoint as it is earlier than the last check point in the queue.\n"
611 << "Checkpoints can only be inserted in ascending order." << std::endl;
612 return;
613 }
614 }
615
616 checkPoints_.push(t);
617 if (this->verbose())
618 std::cout << "Set check point at t = " << t << " seconds." << std::endl;
619 }
620
624 Scalar computeStepSizeRespectingCheckPoints_() const
625 {
626 using std::min;
627 auto maxDt = std::numeric_limits<Scalar>::max();
628
629 if (periodicCheckPoints_)
630 maxDt = min(maxDt, lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time());
631
632 if (!checkPoints_.empty())
633 maxDt = min(maxDt, checkPoints_.front() - this->time());
634
635 return maxDt;
636 }
637
638 bool periodicCheckPoints_;
639 Scalar deltaPeriodicCheckPoint_;
640 Scalar lastPeriodicCheckPoint_;
641 std::queue<Scalar> checkPoints_;
642 bool isCheckPoint_;
643};
644
645} // end namespace Dumux
646
647#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: adapt.hh:29
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
virtual ~TimeLoopBase()
Abstract base class needs virtual constructor.
Definition: timeloop.hh:75
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: timeloop.hh:117
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: timeloop.hh:333
void setMaxTimeStepSize(Scalar maxDt)
Set the maximum time step size to a given value.
Definition: timeloop.hh:271
void resetTimer()
Reset the timer.
Definition: timeloop.hh:150
TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: timeloop.hh:119
void setTime(Scalar t, int stepIdx)
Set the current simulated time and the time step index.
Definition: timeloop.hh:214
Scalar time_
Definition: timeloop.hh:396
int timeStepIndex() const
Returns number of time steps which have been executed since the beginning of the simulation.
Definition: timeloop.hh:289
double stop()
Tells the time loop to stop tracking the time.
Definition: timeloop.hh:142
Scalar endTime_
Definition: timeloop.hh:397
Scalar previousTimeStepSize_
Definition: timeloop.hh:400
void advanceTimeStep() override
Advance time step.
Definition: timeloop.hh:184
int timeStepIdx_
Definition: timeloop.hh:403
void start()
Tells the time loop to start tracking the time.
Definition: timeloop.hh:133
Scalar userSetMaxTimeStepSize_
Definition: timeloop.hh:401
void setFinished(bool finished=true)
Specify whether the simulation is finished.
Definition: timeloop.hh:305
bool verbose() const
If the time loop has verbose output.
Definition: timeloop.hh:387
bool verbose_
Definition: timeloop.hh:405
bool willBeFinished() const
Returns true if the simulation is finished after the time level is incremented by the current time st...
Definition: timeloop.hh:323
Scalar previousTimeStepSize() const
The previous time step size.
Definition: timeloop.hh:295
void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Reset the time loop.
Definition: timeloop.hh:158
Scalar endTime() const
Returns the number of (simulated) seconds which the simulation runs.
Definition: timeloop.hh:228
void setTimeStepSize(Scalar dt) final
Set the current time step size to a given value.
Definition: timeloop.hh:259
bool finished() const override
Returns true if the simulation is finished.
Definition: timeloop.hh:314
Scalar timeStepSize_
Definition: timeloop.hh:399
void finalize(const Communicator &comm=Dune::MPIHelper::getCollectiveCommunication())
Print final status and stops tracking the time.
Definition: timeloop.hh:367
void reportTimeStep() const
State info on cpu time.
Definition: 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: timeloop.hh:282
bool finished_
Definition: timeloop.hh:404
Dune::Timer timer_
Definition: timeloop.hh:395
void setTime(Scalar t)
Set the current simulated time, don't change the current time step index.
Definition: timeloop.hh:205
void setEndTime(Scalar t)
Set the time of simulated seconds at which the simulation runs.
Definition: timeloop.hh:236
double wallClockTime() const
Returns the current wall clock time (cpu time) spend in this time loop.
Definition: timeloop.hh:246
Scalar timeAfterLastTimeStep_
Definition: timeloop.hh:402
Scalar time() const final
Return the time before the time integration. To get the time after the time integration you have to ...
Definition: timeloop.hh:222
Scalar timeStepWallClockTime_
Definition: timeloop.hh:402
A time loop with a check point mechanism.
Definition: timeloop.hh:414
bool isCheckPoint() const
Whether now is a time checkpoint.
Definition: timeloop.hh:548
void setPeriodicCheckPoint(Scalar interval, Scalar offset=0.0)
Set a periodic check point.
Definition: timeloop.hh:508
void setCheckPoint(Scalar t)
add a checkpoint to the queue
Definition: timeloop.hh:557
void advanceTimeStep() override
Advance time step.
Definition: timeloop.hh:428
void removeAllCheckPoints()
remove all check points
Definition: timeloop.hh:537
void setCheckPoint(const std::vector< Scalar > &checkPoints)
add checkpoints to the queue from a vector of time points
Definition: timeloop.hh:572
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: timeloop.hh:490
void setCheckPoint(ForwardIterator first, ForwardIterator last)
add checkpoints to the queue from a container from the first iterator to the last iterator
Definition: timeloop.hh:583
CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: timeloop.hh:416
void disablePeriodicCheckPoints()
disable periodic check points
Definition: timeloop.hh:533