3.1-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#include <dune/common/parallel/collectivecommunication.hh>
34#include <dune/common/parallel/mpihelper.hh>
35#include <dune/common/exceptions.hh>
36
38
39namespace Dumux {
40
63template<class Scalar>
65{
66public:
68 virtual ~TimeLoopBase() {};
69
75 virtual Scalar time() const = 0;
76
80 virtual Scalar timeStepSize() const = 0;
81
85 virtual Scalar maxTimeStepSize() const = 0;
86
90 virtual void advanceTimeStep() = 0;
91
96 virtual void setTimeStepSize(Scalar dt) = 0;
97
101 virtual bool finished() const = 0;
102};
103
108template <class Scalar>
109class TimeLoop : public TimeLoopBase<Scalar>
110{
111public:
112 TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
113 : timer_(false)
114 {
115 reset(startTime, dt, tEnd, verbose);
116 }
117
126 void start()
127 {
128 timer_.start();
129 }
130
135 double stop()
136 {
137 return timer_.stop();
138 }
139
144 {
145 timer_.reset();
146 }
147
151 void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
152 {
153 verbose_ =
154 verbose &&
155 Dune::MPIHelper::getCollectiveCommunication().rank() == 0;
156
157 time_ = startTime;
158 endTime_ = tEnd;
159
161 userSetMaxTimeStepSize_ = std::numeric_limits<Scalar>::max();
162 timeStepIdx_ = 0;
163 finished_ = false;
166
167 // ensure that dt is not greater than tEnd-startTime
168 setTimeStepSize(dt);
169
170 timer_.stop();
171 timer_.reset();
172 }
173
177 void advanceTimeStep() override
178 {
179 timeStepIdx_++;
182
183 // compute how long the last time step took
184 const auto cpuTime = wallClockTime();
186 timeAfterLastTimeStep_ = cpuTime;
187
188 // ensure that using current dt we don't exceed tEnd in next time step
190 }
191
198 void setTime(Scalar t)
199 { time_ = t; }
200
207 void setTime(Scalar t, int stepIdx)
208 { time_ = t; timeStepIdx_ = stepIdx; }
209
215 Scalar time() const final
216 { return time_; }
217
221 Scalar endTime() const
222 { return endTime_; }
223
229 void setEndTime(Scalar t)
230 {
231 endTime_ = t;
232 if (verbose_)
233 std::cout << "Set new end time to t = " << t << " seconds." << std::endl;
234 }
235
239 double wallClockTime() const
240 { return timer_.elapsed(); }
241
252 void setTimeStepSize(Scalar dt) final
253 {
254 using std::min;
255 timeStepSize_ = min(dt, maxTimeStepSize());
256 }
257
264 void setMaxTimeStepSize(Scalar maxDt)
265 {
268 }
269
275 Scalar timeStepSize() const final
276 { return timeStepSize_; }
277
282 int timeStepIndex() const
283 { return timeStepIdx_; }
284
288 Scalar previousTimeStepSize() const
289 { return previousTimeStepSize_; }
290
298 void setFinished(bool finished = true)
299 { finished_ = finished; }
300
307 bool finished() const override
308 {
309 return finished_ || endTime_-time_ < 1e-10*time_;
310 }
311
316 bool willBeFinished() const
317 {
319 }
320
326 Scalar maxTimeStepSize() const override
327 {
328 if (finished())
329 return 0.0;
330
331 using std::min; using std::max;
332 return min(userSetMaxTimeStepSize_, max<Scalar>(0.0, endTime_ - time_));
333 }
334
339 void reportTimeStep() const
340 {
341 if (verbose_)
342 {
343 const auto cpuTime = wallClockTime();
344 const auto percent = std::round( time_ / endTime_ * 100 );
345 std::cout << "[" << std::fixed << std::setw( 3 ) << std::setfill( ' ' )
346 << std::setprecision( 0 ) << percent << "%] "
347 << "Time step " << timeStepIdx_ << " done in "
348 << std::setprecision( 6 ) << timeStepWallClockTime_ << " seconds. "
349 << "Wall clock time: " << std::setprecision( 3 ) << cpuTime
350 << ", time: " << std::setprecision( 5 ) << time_
351 << ", time step size: " << std::setprecision( 8 ) << previousTimeStepSize_
352 << std::endl;
353 }
354 }
355
359 template< class Communicator = Dune::CollectiveCommunication<typename Dune::MPIHelper::MPICommunicator> >
360 void finalize(const Communicator& comm = Dune::MPIHelper::getCollectiveCommunication())
361 {
362 auto cpuTime = timer_.stop();
363
364 if (verbose_)
365 {
366 std::cout << "Simulation took " << cpuTime << " seconds on "
367 << comm.size() << " processes.\n";
368 }
369
370 if (comm.size() > 1)
371 cpuTime = comm.sum(cpuTime);
372
373 if (verbose_)
374 {
375 std::cout << "The cumulative CPU time was " << cpuTime << " seconds.\n";
376 }
377 }
378
380 bool verbose() const
381 { return verbose_; }
382
383 /*
384 * @}
385 */
386
387protected:
388 Dune::Timer timer_;
389 Scalar time_;
390 Scalar endTime_;
391
399};
400
405template <class Scalar>
406class CheckPointTimeLoop : public TimeLoop<Scalar>
407{
408public:
409 CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose = true)
410 : TimeLoop<Scalar>(startTime, dt, tEnd, verbose)
411 {
412 periodicCheckPoints_ = false;
413 deltaPeriodicCheckPoint_ = 0.0;
414 lastPeriodicCheckPoint_ = startTime;
415 isCheckPoint_ = false;
416 }
417
421 void advanceTimeStep() override
422 {
423 const auto dt = this->timeStepSize();
424 const auto newTime = this->time()+dt;
425
427 // if we reached a periodic check point
428 if (periodicCheckPoints_ && Dune::FloatCmp::eq(newTime - lastPeriodicCheckPoint_, deltaPeriodicCheckPoint_, 1e-7))
429 {
430 lastPeriodicCheckPoint_ += deltaPeriodicCheckPoint_;
431 isCheckPoint_ = true;
432 }
433
434 // or a manually set check point
435 else if (!checkPoints_.empty() && Dune::FloatCmp::eq(newTime - checkPoints_.front(), 0.0, 1e-7))
436 {
437 checkPoints_.pop();
438 isCheckPoint_ = true;
439 }
440
441 // if not reset the check point flag
442 else
443 {
444 isCheckPoint_ = false;
445 }
446
447 const auto previousTimeStepSize = this->previousTimeStepSize();
448
449 // advance the time step like in the parent class
451
452 // if this is a check point we might have reduced the time step to reach this check point
453 // reset the time step size to the time step size before this time step
454 if (!this->willBeFinished())
455 {
456 using std::max;
457 if (isCheckPoint_)
459
460 // if there is a check point soon check if the time step after the next time step would be smaller
461 // than 20% of the next time step, if yes increase the suggested next time step to exactly reach the check point
462 // (in the limits of the maximum time step size)
463 auto nextDt = this->timeStepSize();
464 const auto threshold = 0.2*nextDt;
465 const auto nextTime = this->time() + nextDt;
466
467 if (periodicCheckPoints_ && Dune::FloatCmp::le(lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - nextTime, threshold))
468 nextDt = lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time();
469
470 if (!checkPoints_.empty() && Dune::FloatCmp::le(checkPoints_.front() - nextTime, threshold))
471 nextDt = checkPoints_.front() - this->time();
472
473 assert(nextDt > 0.0);
474 this->setTimeStepSize(nextDt);
475 }
476 }
477
483 Scalar maxTimeStepSize() const override
484 {
485 using std::min;
486 const auto maxCheckPointDt = computeStepSizeRespectingCheckPoints_();
487 const auto maxDtParent = TimeLoop<Scalar>::maxTimeStepSize();
488 return min(maxDtParent, maxCheckPointDt);
489 }
490
501 void setPeriodicCheckPoint(Scalar interval, Scalar offset = 0.0)
502 {
503 using std::signbit;
504 if (signbit(interval))
505 DUNE_THROW(Dune::InvalidStateException, "Interval has to be positive!");
506
507 periodicCheckPoints_ = true;
508 deltaPeriodicCheckPoint_ = interval;
509 lastPeriodicCheckPoint_ = offset;
510 while (lastPeriodicCheckPoint_ + interval - this->time() < 1e-14*interval)
511 lastPeriodicCheckPoint_ += interval;
512
513 if (this->verbose())
514 std::cout << "Enabled periodic check points every " << interval
515 << " seconds with the next check point at " << lastPeriodicCheckPoint_ + interval << " seconds." << std::endl;
516
517 // check if the current time point is a check point
518 if (Dune::FloatCmp::eq(this->time()-lastPeriodicCheckPoint_, 0.0, 1e-7))
519 isCheckPoint_ = true;
520
521 // make sure we respect this check point on the next time step
522 this->setTimeStepSize(this->timeStepSize());
523 }
524
527 { periodicCheckPoints_ = false; }
528
531 {
532 periodicCheckPoints_ = false;
533 while (!checkPoints_.empty())
534 checkPoints_.pop();
535 }
536
541 bool isCheckPoint() const
542 { return isCheckPoint_; }
543
550 void setCheckPoint(Scalar t)
551 {
552 // set the check point
553 setCheckPoint_(t);
554
555 // make sure we respect this check point on the next time step
556 this->setTimeStepSize(this->timeStepSize());
557 }
558
565 void setCheckPoint(const std::vector<Scalar>& checkPoints)
566 { setCheckPoint(checkPoints.begin(), checkPoints.end()); }
567
575 template<class ForwardIterator>
576 void setCheckPoint(ForwardIterator first, ForwardIterator last)
577 {
578 // set the check points
579 for (; first != last; ++first)
580 setCheckPoint_(*first);
581
582 // make sure we respect this check point on the next time step
583 this->setTimeStepSize(this->timeStepSize());
584 }
585
586private:
588 void setCheckPoint_(Scalar t)
589 {
590 if (t < this->time())
591 {
592 if (this->verbose())
593 std::cerr << "Couldn't insert checkpoint at t = " << t
594 << " because that's in the past! (current simulation time is " << this->time() << ")" << std::endl;
595 return;
596 }
597
598 if (!checkPoints_.empty())
599 {
600 if (t < checkPoints_.back())
601 {
602 if (this->verbose())
603 std::cerr << "Couldn't insert checkpoint as it is earlier than the last check point in the queue.\n"
604 << "Checkpoints can only be inserted in ascending order." << std::endl;
605 return;
606 }
607 }
608
609 checkPoints_.push(t);
610 if (this->verbose())
611 std::cout << "Set check point at t = " << t << " seconds." << std::endl;
612 }
613
617 Scalar computeStepSizeRespectingCheckPoints_() const
618 {
619 using std::min;
620 auto maxDt = std::numeric_limits<Scalar>::max();
621
622 if (periodicCheckPoints_)
623 maxDt = min(maxDt, lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_ - this->time());
624
625 if (!checkPoints_.empty())
626 maxDt = min(maxDt, checkPoints_.front() - this->time());
627
628 return maxDt;
629 }
630
631 bool periodicCheckPoints_;
632 Scalar deltaPeriodicCheckPoint_;
633 Scalar lastPeriodicCheckPoint_;
634 std::queue<Scalar> checkPoints_;
635 bool isCheckPoint_;
636};
637
638} // end namespace Dumux
639
640#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Manages the handling of time dependent problems.
Definition: timeloop.hh:65
virtual ~TimeLoopBase()
Abstract base class needs virtual constructor.
Definition: timeloop.hh:68
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:110
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: timeloop.hh:326
void setMaxTimeStepSize(Scalar maxDt)
Set the maximum time step size to a given value.
Definition: timeloop.hh:264
void resetTimer()
Reset the timer.
Definition: timeloop.hh:143
TimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: timeloop.hh:112
void setTime(Scalar t, int stepIdx)
Set the current simulated time and the time step index.
Definition: timeloop.hh:207
Scalar time_
Definition: timeloop.hh:389
int timeStepIndex() const
Returns number of time steps which have been executed since the beginning of the simulation.
Definition: timeloop.hh:282
double stop()
Tells the time loop to stop tracking the time.
Definition: timeloop.hh:135
Scalar endTime_
Definition: timeloop.hh:390
Scalar previousTimeStepSize_
Definition: timeloop.hh:393
void advanceTimeStep() override
Advance time step.
Definition: timeloop.hh:177
int timeStepIdx_
Definition: timeloop.hh:396
void start()
Tells the time loop to start tracking the time.
Definition: timeloop.hh:126
Scalar userSetMaxTimeStepSize_
Definition: timeloop.hh:394
void setFinished(bool finished=true)
Specify whether the simulation is finished.
Definition: timeloop.hh:298
bool verbose() const
If the time loop has verbose output.
Definition: timeloop.hh:380
bool verbose_
Definition: timeloop.hh:398
bool willBeFinished() const
Returns true if the simulation is finished after the time level is incremented by the current time st...
Definition: timeloop.hh:316
Scalar previousTimeStepSize() const
The previous time step size.
Definition: timeloop.hh:288
void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Reset the time loop.
Definition: timeloop.hh:151
Scalar endTime() const
Returns the number of (simulated) seconds which the simulation runs.
Definition: timeloop.hh:221
void setTimeStepSize(Scalar dt) final
Set the current time step size to a given value.
Definition: timeloop.hh:252
bool finished() const override
Returns true if the simulation is finished.
Definition: timeloop.hh:307
Scalar timeStepSize_
Definition: timeloop.hh:392
void finalize(const Communicator &comm=Dune::MPIHelper::getCollectiveCommunication())
Print final status and stops tracking the time.
Definition: timeloop.hh:360
void reportTimeStep() const
State info on cpu time.
Definition: timeloop.hh:339
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:275
bool finished_
Definition: timeloop.hh:397
Dune::Timer timer_
Definition: timeloop.hh:388
void setTime(Scalar t)
Set the current simulated time, don't change the current time step index.
Definition: timeloop.hh:198
void setEndTime(Scalar t)
Set the time of simulated seconds at which the simulation runs.
Definition: timeloop.hh:229
double wallClockTime() const
Returns the current wall clock time (cpu time) spend in this time loop.
Definition: timeloop.hh:239
Scalar timeAfterLastTimeStep_
Definition: timeloop.hh:395
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:215
Scalar timeStepWallClockTime_
Definition: timeloop.hh:395
A time loop with a check point mechanism.
Definition: timeloop.hh:407
bool isCheckPoint() const
Whether now is a time checkpoint.
Definition: timeloop.hh:541
void setPeriodicCheckPoint(Scalar interval, Scalar offset=0.0)
Set a periodic check point.
Definition: timeloop.hh:501
void setCheckPoint(Scalar t)
add a checkpoint to the queue
Definition: timeloop.hh:550
void advanceTimeStep() override
Advance time step.
Definition: timeloop.hh:421
void removeAllCheckPoints()
remove all check points
Definition: timeloop.hh:530
void setCheckPoint(const std::vector< Scalar > &checkPoints)
add checkpoints to the queue from a vector of time points
Definition: timeloop.hh:565
Scalar maxTimeStepSize() const override
The current maximum time step size.
Definition: timeloop.hh:483
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:576
CheckPointTimeLoop(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Definition: timeloop.hh:409
void disablePeriodicCheckPoints()
disable periodic check points
Definition: timeloop.hh:526