version 3.7
multistagemethods.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//
13#ifndef DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
14#define DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
15
16#include <cmath>
17#include <string>
18#include <array>
19
20namespace Dumux::Experimental {
21
57template<class Scalar>
59{
60public:
61 virtual bool implicit () const = 0;
62
63 virtual std::size_t numStages () const = 0;
64
66 virtual Scalar temporalWeight (std::size_t i, std::size_t k) const = 0;
67
69 virtual Scalar spatialWeight (std::size_t i, std::size_t k) const = 0;
70
72 virtual Scalar timeStepWeight (std::size_t k) const = 0;
73
74 virtual std::string name () const = 0;
75
76 virtual ~MultiStageMethod() = default;
77};
78
80namespace MultiStage {
81
88template<class Scalar>
89class Theta : public MultiStageMethod<Scalar>
90{
91public:
92 explicit Theta(const Scalar theta)
93 : paramAlpha_{{-1.0, 1.0}}
94 , paramBeta_{{1.0-theta, theta}}
95 , paramD_{{0.0, 1.0}}
96 {}
97
98 bool implicit () const final
99 { return paramBeta_[1] > 0.0; }
100
101 std::size_t numStages () const final
102 { return 1; }
103
104 Scalar temporalWeight (std::size_t, std::size_t k) const final
105 { return paramAlpha_[k]; }
106
107 Scalar spatialWeight (std::size_t, std::size_t k) const final
108 { return paramBeta_[k]; }
109
110 Scalar timeStepWeight (std::size_t k) const final
111 { return paramD_[k]; }
112
113 std::string name () const override
114 { return "theta scheme"; }
115
116private:
117 std::array<Scalar, 2> paramAlpha_;
118 std::array<Scalar, 2> paramBeta_;
119 std::array<Scalar, 2> paramD_;
120};
121
125template<class Scalar>
126class ExplicitEuler final : public Theta<Scalar>
127{
128public:
129 ExplicitEuler() : Theta<Scalar>(0.0) {}
130
131 std::string name () const final
132 { return "explicit Euler"; }
133};
134
138template<class Scalar>
139class ImplicitEuler final : public Theta<Scalar>
140{
141public:
142 ImplicitEuler() : Theta<Scalar>(1.0) {}
143
144 std::string name () const final
145 { return "implicit Euler"; }
146};
147
151template<class Scalar>
153{
154public:
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}}
165 {}
166
167 bool implicit () const final
168 { return false; }
169
170 std::size_t numStages () const final
171 { return 4; }
172
173 Scalar temporalWeight (std::size_t i, std::size_t k) const final
174 { return paramAlpha_[i-1][k]; }
175
176 Scalar spatialWeight (std::size_t i, std::size_t k) const final
177 { return paramBeta_[i-1][k]; }
178
179 Scalar timeStepWeight (std::size_t k) const final
180 { return paramD_[k]; }
181
182 std::string name () const final
183 { return "explicit Runge-Kutta 4th order"; }
184
185private:
186 std::array<std::array<Scalar, 5>, 4> paramAlpha_;
187 std::array<std::array<Scalar, 5>, 4> paramBeta_;
188 std::array<Scalar, 5> paramD_;
189};
190
191} // end namespace MultiStage
192} // end namespace Dumux::Experimental
193
194#endif
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 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