version 3.7
multistagetimestepper.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_TIMESTEPPING_MULTISTAGE_TIMESTEPPER_HH
13#define DUMUX_TIMESTEPPING_MULTISTAGE_TIMESTEPPER_HH
14
15#include <memory>
16#include <vector>
17#include <cmath>
18
19namespace Dumux::Experimental {
20
22template<class Scalar>
23class MultiStageMethod;
24
26template<class Scalar>
28{
29 struct Params {
30 Scalar alpha, betaDt, timeAtStage, dtFraction;
32 };
33public:
35 MultiStageParams(const MultiStageMethod<Scalar>& m, std::size_t i, const Scalar t, const Scalar dt)
36 : size_(i+1)
37 {
38 params_.resize(size_);
39 for (std::size_t k = 0; k < size_; ++k)
40 {
41 auto& p = params_[k];
42 p.alpha = m.temporalWeight(i, k);
43 p.betaDt = m.spatialWeight(i, k)*dt;
44 p.timeAtStage = t + m.timeStepWeight(k)*dt;
45 p.dtFraction = m.timeStepWeight(k);
46
47 using std::abs;
48 p.skipTemporal = (abs(p.alpha) < 1e-6);
49 p.skipSpatial = (abs(p.betaDt) < 1e-6);
50 }
51 }
52
53 std::size_t size () const
54 { return size_; }
55
57 Scalar temporalWeight (std::size_t k) const
58 { return params_[k].alpha; }
59
61 Scalar spatialWeight (std::size_t k) const
62 { return params_[k].betaDt; }
63
65 Scalar timeAtStage (std::size_t k) const
66 { return params_[k].timeAtStage; }
67
69 Scalar timeStepFraction (std::size_t k) const
70 { return params_[k].dtFraction; }
71
73 Scalar skipTemporal (std::size_t k) const
74 { return params_[k].skipTemporal; }
75
77 Scalar skipSpatial (std::size_t k) const
78 { return params_[k].skipSpatial; }
79
80private:
81 std::size_t size_;
82 std::vector<Params> params_;
83};
84
91template<class PDESolver>
93{
94 using Variables = typename PDESolver::Variables;
95 using Scalar = typename Variables::Scalar;
97
98public:
99
106 MultiStageTimeStepper(std::shared_ptr<PDESolver> pdeSolver,
107 std::shared_ptr<const MultiStageMethod<Scalar>> msMethod)
108 : pdeSolver_(pdeSolver)
109 , msMethod_(msMethod)
110 {}
111
120 void step(Variables& vars, const Scalar t, const Scalar dt)
121 {
122 // make sure there are no traces of previous stages
123 pdeSolver_->assembler().clearStages();
124
125 for (auto stageIdx = 1UL; stageIdx <= msMethod_->numStages(); ++stageIdx)
126 {
127 // extract parameters for this stage from the time stepping method
128 auto stageParams = std::make_shared<StageParams>(*msMethod_, stageIdx, t, dt);
129
130 // prepare the assembler for this stage
131 pdeSolver_->assembler().prepareStage(vars, stageParams);
132
133 // assemble & solve
134 pdeSolver_->solve(vars);
135 }
136
137 // clear traces of previously registered stages
138 pdeSolver_->assembler().clearStages();
139 }
140
144 void setMethod(std::shared_ptr<const MultiStageMethod<Scalar>> msMethod)
145 { msMethod_ = msMethod; }
146
147private:
148 std::shared_ptr<PDESolver> pdeSolver_;
149 std::shared_ptr<const MultiStageMethod<Scalar>> msMethod_;
150};
151
152} // end namespace Dumux::Experimental
153
154#endif
Abstract interface for one-step multi-stage method parameters in Shu/Osher form.
Definition: multistagemethods.hh:59
virtual Scalar timeStepWeight(std::size_t k) const =0
time step weights for each stage ( )
virtual Scalar temporalWeight(std::size_t i, std::size_t k) const =0
weights of the temporal operator residual ( )
virtual Scalar spatialWeight(std::size_t i, std::size_t k) const =0
weights of the spatial operator residual ( )
Data object for the parameters of a given stage.
Definition: multistagetimestepper.hh:28
Scalar timeAtStage(std::size_t k) const
the time at which we have to evaluate the operators
Definition: multistagetimestepper.hh:65
Scalar spatialWeight(std::size_t k) const
weights of the spatial operator residual ( )
Definition: multistagetimestepper.hh:61
Scalar timeStepFraction(std::size_t k) const
the fraction of a time step corresponding to the k-th stage
Definition: multistagetimestepper.hh:69
Scalar skipTemporal(std::size_t k) const
If .
Definition: multistagetimestepper.hh:73
std::size_t size() const
Definition: multistagetimestepper.hh:53
Scalar skipSpatial(std::size_t k) const
If .
Definition: multistagetimestepper.hh:77
Scalar temporalWeight(std::size_t k) const
weights of the temporal operator residual ( )
Definition: multistagetimestepper.hh:57
MultiStageParams(const MultiStageMethod< Scalar > &m, std::size_t i, const Scalar t, const Scalar dt)
Extract params for stage i from method m.
Definition: multistagetimestepper.hh:35
Time stepping with a multi-stage method.
Definition: multistagetimestepper.hh:93
void setMethod(std::shared_ptr< const MultiStageMethod< Scalar > > msMethod)
Set/change the time step method.
Definition: multistagetimestepper.hh:144
void step(Variables &vars, const Scalar t, const Scalar dt)
Advance one time step of the given time loop.
Definition: multistagetimestepper.hh:120
MultiStageTimeStepper(std::shared_ptr< PDESolver > pdeSolver, std::shared_ptr< const MultiStageMethod< Scalar > > msMethod)
The constructor.
Definition: multistagetimestepper.hh:106
typename ScalarT< X >::type Scalar
export the underlying scalar type
Definition: variables.hh:57
Detail::PDESolver::AssemblerVariables< Assembler > Variables
export the type of variables that represent a numerical solution
Definition: common/pdesolver.hh:71
Definition: variables.hh:21