version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © dune-pdelab developers, see LICENSE.md (permalink below)
7// SPDX-License-Identifier: GPL-3.0-or-later
8// The code is based on the implementation of time stepping method parameters
9// in dune-pdelab (https://archive.softwareheritage.org/swh:1:cnt:9c3d72412f8e4d48d84a0090b2bb4362b5a0d843)
10// licensed under GPL-3.0-or-later, see their LICENSE.md for a full list of copyright holders at
11// https://archive.softwareheritage.org/swh:1:cnt:b11b484e74eefe20c74f7043309b1b02853df2eb.
12// Modifications (different interface naming, comments, types) are licensed under GPL-3.0-or-later.
13//
20#ifndef DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
21#define DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
22
23#include <cmath>
24#include <string>
25#include <array>
26
27namespace Dumux::Experimental {
28
73template<class Scalar>
75{
76public:
77 virtual bool implicit () const = 0;
78
79 virtual std::size_t numStages () const = 0;
80
82 virtual Scalar temporalWeight (std::size_t i, std::size_t k) const = 0;
83
85 virtual Scalar spatialWeight (std::size_t i, std::size_t k) const = 0;
86
88 virtual Scalar timeStepWeight (std::size_t k) const = 0;
89
90 virtual std::string name () const = 0;
91
92 virtual ~MultiStageMethod() = default;
93};
94
96namespace MultiStage {
97
104template<class Scalar>
105class Theta : public MultiStageMethod<Scalar>
106{
107public:
108 explicit Theta(const Scalar theta)
109 : paramAlpha_{{-1.0, 1.0}}
110 , paramBeta_{{1.0-theta, theta}}
111 , paramD_{{0.0, 1.0}}
112 {}
113
114 bool implicit () const final
115 { return paramBeta_[1] > 0.0; }
116
117 std::size_t numStages () const final
118 { return 1; }
119
120 Scalar temporalWeight (std::size_t, std::size_t k) const final
121 { return paramAlpha_[k]; }
122
123 Scalar spatialWeight (std::size_t, std::size_t k) const final
124 { return paramBeta_[k]; }
125
126 Scalar timeStepWeight (std::size_t k) const final
127 { return paramD_[k]; }
128
129 std::string name () const override
130 { return "theta scheme"; }
131
132private:
133 std::array<Scalar, 2> paramAlpha_;
134 std::array<Scalar, 2> paramBeta_;
135 std::array<Scalar, 2> paramD_;
136};
137
141template<class Scalar>
142class ExplicitEuler final : public Theta<Scalar>
143{
144public:
145 ExplicitEuler() : Theta<Scalar>(0.0) {}
146
147 std::string name () const final
148 { return "explicit Euler"; }
149};
150
154template<class Scalar>
155class ImplicitEuler final : public Theta<Scalar>
156{
157public:
158 ImplicitEuler() : Theta<Scalar>(1.0) {}
159
160 std::string name () const final
161 { return "implicit Euler"; }
162};
163
167template<class Scalar>
169{
170public:
172 : paramAlpha_{{{-1.0, 1.0, 0.0, 0.0, 0.0},
173 {-1.0, 0.0, 1.0, 0.0, 0.0},
174 {-1.0, 0.0, 0.0, 1.0, 0.0},
175 {-1.0, 0.0, 0.0, 0.0, 1.0}}}
176 , paramBeta_{{{0.5, 0.0, 0.0, 0.0, 0.0},
177 {0.0, 0.5, 0.0, 0.0, 0.0},
178 {0.0, 0.0, 1.0, 0.0, 0.0},
179 {1.0/6.0, 1.0/3.0, 1.0/3.0, 1.0/6.0, 0.0}}}
180 , paramD_{{0.0, 0.5, 0.5, 1.0, 1.0}}
181 {}
182
183 bool implicit () const final
184 { return false; }
185
186 std::size_t numStages () const final
187 { return 4; }
188
189 Scalar temporalWeight (std::size_t i, std::size_t k) const final
190 { return paramAlpha_[i-1][k]; }
191
192 Scalar spatialWeight (std::size_t i, std::size_t k) const final
193 { return paramBeta_[i-1][k]; }
194
195 Scalar timeStepWeight (std::size_t k) const final
196 { return paramD_[k]; }
197
198 std::string name () const final
199 { return "explicit Runge-Kutta 4th order"; }
200
201private:
202 std::array<std::array<Scalar, 5>, 4> paramAlpha_;
203 std::array<std::array<Scalar, 5>, 4> paramBeta_;
204 std::array<Scalar, 5> paramD_;
205};
206
211template<class Scalar>
212class DIRKThirdOrderAlexander final : public MultiStageMethod<Scalar>
213{
214public:
216 {
217 constexpr Scalar alpha = []{
218 // Newton iteration for alpha
219 Scalar alpha = 0.4358665215; // initial guess
220 for (int i = 0; i < 10; ++i)
221 alpha = alpha - (alpha*(alpha*alpha-3.0*(alpha-0.5))-1.0/6.0)/(3.0*alpha*(alpha-2.0)+1.5);
222 return alpha;
223 }();
224
225 constexpr Scalar tau2 = (1.0+alpha)*0.5;
226 constexpr Scalar b1 = -(6.0*alpha*alpha -16.0*alpha + 1.0)*0.25;
227 constexpr Scalar b2 = (6*alpha*alpha - 20.0*alpha + 5.0)*0.25;
228
229 paramD_ = {{0.0, alpha, tau2, 1.0}};
230 paramAlpha_ = {{
231 {-1.0, 1.0, 0.0, 0.0},
232 {-1.0, 0.0, 1.0, 0.0},
233 {-1.0, 0.0, 0.0, 1.0}
234 }};
235 paramBeta_ = {{
236 {0.0, alpha, 0.0, 0.0},
237 {0.0, tau2-alpha, alpha, 0.0},
238 {0.0, b1, b2, alpha}
239 }};
240 }
241
242 bool implicit () const final
243 { return true; }
244
245 std::size_t numStages () const final
246 { return 3; }
247
248 Scalar temporalWeight (std::size_t i, std::size_t k) const final
249 { return paramAlpha_[i-1][k]; }
250
251 Scalar spatialWeight (std::size_t i, std::size_t k) const final
252 { return paramBeta_[i-1][k]; }
253
254 Scalar timeStepWeight (std::size_t k) const final
255 { return paramD_[k]; }
256
257 std::string name () const final
258 { return "diagonally implicit Runge-Kutta 3rd order (Alexander)"; }
259
260private:
261 std::array<std::array<Scalar, 4>, 3> paramAlpha_;
262 std::array<std::array<Scalar, 4>, 3> paramBeta_;
263 std::array<Scalar, 4> paramD_;
264};
265
266} // end namespace MultiStage
267} // end namespace Dumux::Experimental
268
269#endif
Third order DIRK scheme.
Definition: multistagemethods.hh:213
bool implicit() const final
Definition: multistagemethods.hh:242
DIRKThirdOrderAlexander()
Definition: multistagemethods.hh:215
Scalar temporalWeight(std::size_t i, std::size_t k) const final
weights of the temporal operator residual ( )
Definition: multistagemethods.hh:248
Scalar spatialWeight(std::size_t i, std::size_t k) const final
weights of the spatial operator residual ( )
Definition: multistagemethods.hh:251
Scalar timeStepWeight(std::size_t k) const final
time step weights for each stage ( )
Definition: multistagemethods.hh:254
std::string name() const final
Definition: multistagemethods.hh:257
std::size_t numStages() const final
Definition: multistagemethods.hh:245
An explicit Euler time stepping scheme.
Definition: multistagemethods.hh:143
ExplicitEuler()
Definition: multistagemethods.hh:145
std::string name() const final
Definition: multistagemethods.hh:147
An implicit Euler time stepping scheme.
Definition: multistagemethods.hh:156
ImplicitEuler()
Definition: multistagemethods.hh:158
std::string name() const final
Definition: multistagemethods.hh:160
Classical explicit fourth order Runge-Kutta scheme.
Definition: multistagemethods.hh:169
RungeKuttaExplicitFourthOrder()
Definition: multistagemethods.hh:171
Scalar temporalWeight(std::size_t i, std::size_t k) const final
weights of the temporal operator residual ( )
Definition: multistagemethods.hh:189
std::string name() const final
Definition: multistagemethods.hh:198
bool implicit() const final
Definition: multistagemethods.hh:183
std::size_t numStages() const final
Definition: multistagemethods.hh:186
Scalar spatialWeight(std::size_t i, std::size_t k) const final
weights of the spatial operator residual ( )
Definition: multistagemethods.hh:192
Scalar timeStepWeight(std::size_t k) const final
time step weights for each stage ( )
Definition: multistagemethods.hh:195
A theta time stepping scheme theta=1.0 is an implicit Euler scheme, theta=0.0 an explicit Euler schem...
Definition: multistagemethods.hh:106
std::size_t numStages() const final
Definition: multistagemethods.hh:117
bool implicit() const final
Definition: multistagemethods.hh:114
std::string name() const override
Definition: multistagemethods.hh:129
Theta(const Scalar theta)
Definition: multistagemethods.hh:108
Scalar spatialWeight(std::size_t, std::size_t k) const final
weights of the spatial operator residual ( )
Definition: multistagemethods.hh:123
Scalar temporalWeight(std::size_t, std::size_t k) const final
weights of the temporal operator residual ( )
Definition: multistagemethods.hh:120
Scalar timeStepWeight(std::size_t k) const final
time step weights for each stage ( )
Definition: multistagemethods.hh:126
Abstract interface for one-step multi-stage method parameters in Shu/Osher form.
Definition: multistagemethods.hh:75
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: experimental/assembly/cclocalassembler.hh:36