13#ifndef DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
14#define DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
74 virtual std::string
name ()
const = 0;
92 explicit Theta(
const Scalar theta)
93 : paramAlpha_{{-1.0, 1.0}}
94 , paramBeta_{{1.0-theta, theta}}
99 {
return paramBeta_[1] > 0.0; }
105 {
return paramAlpha_[k]; }
108 {
return paramBeta_[k]; }
111 {
return paramD_[k]; }
113 std::string
name ()
const override
114 {
return "theta scheme"; }
117 std::array<Scalar, 2> paramAlpha_;
118 std::array<Scalar, 2> paramBeta_;
119 std::array<Scalar, 2> paramD_;
125template<
class Scalar>
131 std::string
name () const final
132 {
return "explicit Euler"; }
138template<
class Scalar>
144 std::string
name () const final
145 {
return "implicit Euler"; }
151template<
class Scalar>
156 : paramAlpha_{{{-1.0, 1.0, 0.0, 0.0, 0.0},
157 {-1.0, 0.0, 1.0, 0.0, 0.0},
158 {-1.0, 0.0, 0.0, 1.0, 0.0},
159 {-1.0, 0.0, 0.0, 0.0, 1.0}}}
160 , paramBeta_{{{0.5, 0.0, 0.0, 0.0, 0.0},
161 {0.0, 0.5, 0.0, 0.0, 0.0},
162 {0.0, 0.0, 1.0, 0.0, 0.0},
163 {1.0/6.0, 1.0/3.0, 1.0/3.0, 1.0/6.0, 0.0}}}
164 , paramD_{{0.0, 0.5, 0.5, 1.0, 1.0}}
174 {
return paramAlpha_[i-1][k]; }
177 {
return paramBeta_[i-1][k]; }
180 {
return paramD_[k]; }
182 std::string
name () const final
183 {
return "explicit Runge-Kutta 4th order"; }
186 std::array<std::array<Scalar, 5>, 4> paramAlpha_;
187 std::array<std::array<Scalar, 5>, 4> paramBeta_;
188 std::array<Scalar, 5> paramD_;
An explicit Euler time stepping scheme.
Definition: multistagemethods.hh:127
ExplicitEuler()
Definition: multistagemethods.hh:129
std::string name() const final
Definition: multistagemethods.hh:131
An implicit Euler time stepping scheme.
Definition: multistagemethods.hh:140
ImplicitEuler()
Definition: multistagemethods.hh:142
std::string name() const final
Definition: multistagemethods.hh:144
Classical explicit fourth order Runge-Kutta scheme.
Definition: multistagemethods.hh:153
RungeKuttaExplicitFourthOrder()
Definition: multistagemethods.hh:155
Scalar temporalWeight(std::size_t i, std::size_t k) const final
weights of the temporal operator residual ( )
Definition: multistagemethods.hh:173
std::string name() const final
Definition: multistagemethods.hh:182
bool implicit() const final
Definition: multistagemethods.hh:167
std::size_t numStages() const final
Definition: multistagemethods.hh:170
Scalar spatialWeight(std::size_t i, std::size_t k) const final
weights of the spatial operator residual ( )
Definition: multistagemethods.hh:176
Scalar timeStepWeight(std::size_t k) const final
time step weights for each stage ( )
Definition: multistagemethods.hh:179
A theta time stepping scheme theta=1.0 is an implicit Euler scheme, theta=0.0 an explicit Euler schem...
Definition: multistagemethods.hh:90
std::size_t numStages() const final
Definition: multistagemethods.hh:101
bool implicit() const final
Definition: multistagemethods.hh:98
std::string name() const override
Definition: multistagemethods.hh:113
Theta(const Scalar theta)
Definition: multistagemethods.hh:92
Scalar spatialWeight(std::size_t, std::size_t k) const final
weights of the spatial operator residual ( )
Definition: multistagemethods.hh:107
Scalar temporalWeight(std::size_t, std::size_t k) const final
weights of the temporal operator residual ( )
Definition: multistagemethods.hh:104
Scalar timeStepWeight(std::size_t k) const final
time step weights for each stage ( )
Definition: multistagemethods.hh:110
Abstract interface for one-step multi-stage method parameters in Shu/Osher form.
Definition: multistagemethods.hh:59
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 ( )
Definition: variables.hh:21