3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
python/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_PYTHON_COMMON_TIMELOOP_HH
25#define DUMUX_PYTHON_COMMON_TIMELOOP_HH
26
28
29#include <dune/python/pybind11/pybind11.h>
30#include <dune/python/pybind11/stl.h>
31
32namespace Dumux::Python {
33
34template <class Scalar, class... options>
35void registerTimeLoop(pybind11::handle scope,
36 pybind11::class_<CheckPointTimeLoop<Scalar>, options...> cls)
37{
38 using pybind11::operator""_a;
39
41 cls.def(pybind11::init([](Scalar startTime, Scalar dt, Scalar endTime, bool verbose){
42 return new TimeLoop(startTime, dt, endTime, verbose);
43 }), "startTime"_a, "dt"_a, "endTime"_a, "verbose"_a=true);
44
45 cls.def_property_readonly("time", &TimeLoop::time);
46 cls.def_property_readonly("timeStepSize", &TimeLoop::timeStepSize);
47 cls.def_property_readonly("finished", &TimeLoop::finished);
48 cls.def_property_readonly("isCheckPoint", &TimeLoop::isCheckPoint);
49 cls.def("start", &TimeLoop::start);
50 cls.def("reset", &TimeLoop::reset, "startTime"_a, "dt"_a, "endTime"_a, "verbose"_a=true);
51 cls.def("advanceTimeStep", &TimeLoop::advanceTimeStep);
52 cls.def("setTimeStepSize", &TimeLoop::setTimeStepSize, "dt"_a);
53 cls.def("setMaxTimeStepSize", &TimeLoop::setMaxTimeStepSize, "dt"_a);
54 cls.def("reportTimeStep", &TimeLoop::reportTimeStep);
55 cls.def("finalize", [](TimeLoop& self){ self.finalize(); });
56 cls.def("setPeriodicCheckPoint", &TimeLoop::setPeriodicCheckPoint, "interval"_a, "offset"_a=0.0);
57 cls.def("setCheckPoints", [](TimeLoop& self, const std::vector<double>& checkPoints) {
58 self.setCheckPoint(checkPoints);
59 });
60}
61
62template<class Scalar>
63void registerTimeLoop(pybind11::handle scope, const char *clsName = "TimeLoop")
64{
65 pybind11::class_<CheckPointTimeLoop<Scalar>> cls(scope, clsName);
66 registerTimeLoop(scope, cls);
67}
68
69} // namespace Dumux::Python
70
71#endif
Definition: python/assembly/fvassembler.hh:30
void registerTimeLoop(pybind11::handle scope, pybind11::class_< CheckPointTimeLoop< Scalar >, options... > cls)
Definition: python/common/timeloop.hh:35
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
void setMaxTimeStepSize(Scalar maxDt)
Set the maximum time step size to a given value.
Definition: common/timeloop.hh:271
void advanceTimeStep() override
Advance time step.
Definition: common/timeloop.hh:181
void start()
Tells the time loop to start tracking the time.
Definition: common/timeloop.hh:130
void finalize(const Communicator &comm=Dune::MPIHelper::getCommunication())
Print final status and stops tracking the time.
Definition: common/timeloop.hh:363
void reset(Scalar startTime, Scalar dt, Scalar tEnd, bool verbose=true)
Reset the time loop.
Definition: common/timeloop.hh:155
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
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
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
A time loop with a check point mechanism.
Definition: common/timeloop.hh:409
Manages the handling of time dependent problems.