25#ifndef DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
26#define DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
86 virtual std::string
name ()
const = 0;
100template<
class Scalar>
105 : paramAlpha_{{-1.0, 1.0}}
106 , paramBeta_{{1.0-theta, theta}}
107 , paramD_{{0.0, 1.0}}
111 {
return paramBeta_[1] > 0.0; }
117 {
return paramAlpha_[k]; }
120 {
return paramBeta_[k]; }
123 {
return paramD_[k]; }
125 std::string
name ()
const override
126 {
return "theta scheme"; }
129 std::array<Scalar, 2> paramAlpha_;
130 std::array<Scalar, 2> paramBeta_;
131 std::array<Scalar, 2> paramD_;
137template<
class Scalar>
143 std::string
name () const final
144 {
return "explicit Euler"; }
150template<
class Scalar>
156 std::string
name () const final
157 {
return "implicit Euler"; }
163template<
class Scalar>
168 : paramAlpha_{{{-1.0, 1.0, 0.0, 0.0, 0.0},
169 {-1.0, 0.0, 1.0, 0.0, 0.0},
170 {-1.0, 0.0, 0.0, 1.0, 0.0},
171 {-1.0, 0.0, 0.0, 0.0, 1.0}}}
172 , paramBeta_{{{0.5, 0.0, 0.0, 0.0, 0.0},
173 {0.0, 0.5, 0.0, 0.0, 0.0},
174 {0.0, 0.0, 1.0, 0.0, 0.0},
175 {1.0/6.0, 1.0/3.0, 1.0/3.0, 1.0/6.0, 0.0}}}
176 , paramD_{{0.0, 0.5, 0.5, 1.0, 1.0}}
186 {
return paramAlpha_[i-1][k]; }
189 {
return paramBeta_[i-1][k]; }
192 {
return paramD_[k]; }
194 std::string
name () const final
195 {
return "explicit Runge-Kutta 4th order"; }
198 std::array<std::array<Scalar, 5>, 4> paramAlpha_;
199 std::array<std::array<Scalar, 5>, 4> paramBeta_;
200 std::array<Scalar, 5> paramD_;
Definition: variables.hh:32
Abstract interface for one-step multi-stage method parameters in Shu/Osher form.
Definition: multistagemethods.hh:71
virtual std::string name() const =0
virtual std::size_t numStages() const =0
virtual ~MultiStageMethod()=default
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 bool implicit() const =0
virtual Scalar spatialWeight(std::size_t i, std::size_t k) const =0
weights of the spatial operator residual ( )
A theta time stepping scheme theta=1.0 is an implicit Euler scheme, theta=0.0 an explicit Euler schem...
Definition: multistagemethods.hh:102
std::size_t numStages() const final
Definition: multistagemethods.hh:113
bool implicit() const final
Definition: multistagemethods.hh:110
std::string name() const override
Definition: multistagemethods.hh:125
Theta(const Scalar theta)
Definition: multistagemethods.hh:104
Scalar spatialWeight(std::size_t, std::size_t k) const final
weights of the spatial operator residual ( )
Definition: multistagemethods.hh:119
Scalar temporalWeight(std::size_t, std::size_t k) const final
weights of the temporal operator residual ( )
Definition: multistagemethods.hh:116
Scalar timeStepWeight(std::size_t k) const final
time step weights for each stage ( )
Definition: multistagemethods.hh:122
An explicit Euler time stepping scheme.
Definition: multistagemethods.hh:139
ExplicitEuler()
Definition: multistagemethods.hh:141
std::string name() const final
Definition: multistagemethods.hh:143
An implicit Euler time stepping scheme.
Definition: multistagemethods.hh:152
ImplicitEuler()
Definition: multistagemethods.hh:154
std::string name() const final
Definition: multistagemethods.hh:156
Classical explicit fourth order Runge-Kutta scheme.
Definition: multistagemethods.hh:165
RungeKuttaExplicitFourthOrder()
Definition: multistagemethods.hh:167
Scalar temporalWeight(std::size_t i, std::size_t k) const final
weights of the temporal operator residual ( )
Definition: multistagemethods.hh:185
std::string name() const final
Definition: multistagemethods.hh:194
bool implicit() const final
Definition: multistagemethods.hh:179
std::size_t numStages() const final
Definition: multistagemethods.hh:182
Scalar spatialWeight(std::size_t i, std::size_t k) const final
weights of the spatial operator residual ( )
Definition: multistagemethods.hh:188
Scalar timeStepWeight(std::size_t k) const final
time step weights for each stage ( )
Definition: multistagemethods.hh:191