3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
26#define DUMUX_TIMESTEPPING_MULTISTAGE_METHODS_HH
27
28#include <cmath>
29#include <string>
30#include <array>
31
32namespace Dumux::Experimental {
33
69template<class Scalar>
71{
72public:
73 virtual bool implicit () const = 0;
74
75 virtual std::size_t numStages () const = 0;
76
78 virtual Scalar temporalWeight (std::size_t i, std::size_t k) const = 0;
79
81 virtual Scalar spatialWeight (std::size_t i, std::size_t k) const = 0;
82
84 virtual Scalar timeStepWeight (std::size_t k) const = 0;
85
86 virtual std::string name () const = 0;
87
88 virtual ~MultiStageMethod() = default;
89};
90
92namespace MultiStage {
93
100template<class Scalar>
101class Theta : public MultiStageMethod<Scalar>
102{
103public:
104 explicit Theta(const Scalar theta)
105 : paramAlpha_{{-1.0, 1.0}}
106 , paramBeta_{{1.0-theta, theta}}
107 , paramD_{{0.0, 1.0}}
108 {}
109
110 bool implicit () const final
111 { return paramBeta_[1] > 0.0; }
112
113 std::size_t numStages () const final
114 { return 1; }
115
116 Scalar temporalWeight (std::size_t, std::size_t k) const final
117 { return paramAlpha_[k]; }
118
119 Scalar spatialWeight (std::size_t, std::size_t k) const final
120 { return paramBeta_[k]; }
121
122 Scalar timeStepWeight (std::size_t k) const final
123 { return paramD_[k]; }
124
125 std::string name () const override
126 { return "theta scheme"; }
127
128private:
129 std::array<Scalar, 2> paramAlpha_;
130 std::array<Scalar, 2> paramBeta_;
131 std::array<Scalar, 2> paramD_;
132};
133
137template<class Scalar>
138class ExplicitEuler final : public Theta<Scalar>
139{
140public:
141 ExplicitEuler() : Theta<Scalar>(0.0) {}
142
143 std::string name () const final
144 { return "explicit Euler"; }
145};
146
150template<class Scalar>
151class ImplicitEuler final : public Theta<Scalar>
152{
153public:
154 ImplicitEuler() : Theta<Scalar>(1.0) {}
155
156 std::string name () const final
157 { return "implicit Euler"; }
158};
159
163template<class Scalar>
165{
166public:
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}}
177 {}
178
179 bool implicit () const final
180 { return false; }
181
182 std::size_t numStages () const final
183 { return 4; }
184
185 Scalar temporalWeight (std::size_t i, std::size_t k) const final
186 { return paramAlpha_[i-1][k]; }
187
188 Scalar spatialWeight (std::size_t i, std::size_t k) const final
189 { return paramBeta_[i-1][k]; }
190
191 Scalar timeStepWeight (std::size_t k) const final
192 { return paramD_[k]; }
193
194 std::string name () const final
195 { return "explicit Runge-Kutta 4th order"; }
196
197private:
198 std::array<std::array<Scalar, 5>, 4> paramAlpha_;
199 std::array<std::array<Scalar, 5>, 4> paramBeta_;
200 std::array<Scalar, 5> paramD_;
201};
202
203} // end namespace MultiStage
204} // end namespace Dumux::Experimental
205
206#endif
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 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