version 3.8
staggeredupwindmethods.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_UPWINDING_METHODS_HH
14#define DUMUX_UPWINDING_METHODS_HH
15
16#include <cmath>
17#include <functional>
18#include <iostream>
19#include <string>
20
23
24namespace Dumux {
25
30enum class TvdApproach
31{
32 none, uniform, li, hou
33};
34
40{
42};
43
48template<class Scalar, int upwindSchemeOrder>
50{
51public:
52 StaggeredUpwindMethods(const std::string& paramGroup = "")
53 {
54 upwindWeight_ = getParamFromGroup<Scalar>(paramGroup, "Flux.UpwindWeight");
55
56 if (upwindSchemeOrder > 1)
57 {
58 // Read the runtime parameters
59 tvdApproach_ = tvdApproachFromString(getParamFromGroup<std::string>(paramGroup, "Flux.TvdApproach", "Uniform"));
60 differencingScheme_ = differencingSchemeFromString(getParamFromGroup<std::string>(paramGroup, "Flux.DifferencingScheme", "Minmod"));
61
62 // Assign the limiter_ depending on the differencing scheme
63 switch (differencingScheme_)
64 {
66 {
67 limiter_ = this->vanleer;
68 break;
69 }
71 {
72 if (tvdApproach_ == TvdApproach::hou)
73 DUNE_THROW(ParameterException, "\nDifferencing scheme (Van Albada) is not implemented for the Tvd approach (Hou).");
74 else
75 limiter_ = this->vanalbada;
76 break;
77 }
79 {
80 limiter_ = this->minmod;
81 break;
82 }
84 {
85 limiter_ = this->superbee;
86 break;
87 }
89 {
90 limiter_ = this->umist;
91 break;
92 }
94 {
95 limiter_ = this->mclimiter;
96 break;
97 }
99 {
100 limiter_ = this->wahyd;
101 break;
102 }
103 default:
104 {
105 DUNE_THROW(ParameterException, "\nDifferencing scheme " << static_cast<std::string>(differencingSchemeToString(differencingScheme_)) <<
106 " is not implemented.\n"); // should never be reached
107 break;
108 }
109 }
110
111 if (!hasParamInGroup(paramGroup, "Flux.TvdApproach"))
112 {
113 std::cout << "No TvdApproach specified. Defaulting to the Uniform method." << "\n";
114 std::cout << "Other available TVD approaches for uniform (and nonuniform) grids are as follows: \n"
115 << " " << tvdApproachToString(TvdApproach::uniform) << ": assumes a Uniform cell size distribution\n"
116 << " " << tvdApproachToString(TvdApproach::li) << ": Li's approach for nonuniform cell sizes\n"
117 << " " << tvdApproachToString(TvdApproach::hou) << ": Hou's approach for nonuniform cell sizes \n";
118 std::cout << "Each approach can be specified as written above in the Flux group under the title TvdApproach in your input file. \n";
119 }
120 if (!hasParamInGroup(paramGroup, "Flux.DifferencingScheme"))
121 {
122
123 std::cout << "No DifferencingScheme specified. Defaulting to the Minmod scheme." << "\n";
124 std::cout << "Other available Differencing Schemes are as follows: \n"
125 << " " << differencingSchemeToString(DifferencingScheme::vanleer) << ": The Vanleer flux limiter\n"
126 << " " << differencingSchemeToString(DifferencingScheme::vanalbada) << ": The Vanalbada flux limiter\n"
127 << " " << differencingSchemeToString(DifferencingScheme::minmod) << ": The Minmod flux limiter\n"
128 << " " << differencingSchemeToString(DifferencingScheme::superbee) << ": The Superbee flux limiter\n"
129 << " " << differencingSchemeToString(DifferencingScheme::umist) << ": The Umist flux limiter\n"
130 << " " << differencingSchemeToString(DifferencingScheme::mclimiter) << ": The Mclimiter flux limiter\n"
131 << " " << differencingSchemeToString(DifferencingScheme::wahyd) << ": The Wahyd flux limiter";
132 std::cout << "Each scheme can be specified as written above in the Flux group under the variable DifferencingScheme in your input file. \n";
133 }
134
135 std::cout << "Using the tvdApproach \"" << tvdApproachToString(tvdApproach_)
136 << "\" and the differencing Scheme \" " << differencingSchemeToString(differencingScheme_) << "\" \n";
137 }
138 else
139 {
140 // If the runtime parameters are not specified we will use upwind
141 tvdApproach_ = TvdApproach::none;
142 differencingScheme_ = DifferencingScheme::none;
143 }
144 }
145
151 {
152 if (tvd == "Uniform") return TvdApproach::uniform;
153 if (tvd == "Li") return TvdApproach::li;
154 if (tvd == "Hou") return TvdApproach::hou;
155 DUNE_THROW(ParameterException, "\nThis tvd approach : \"" << tvd << "\" is not implemented.\n"
156 << "The available TVD approaches for uniform (and nonuniform) grids are as follows: \n"
157 << tvdApproachToString(TvdApproach::uniform) << ": assumes a uniform cell size distribution\n"
158 << tvdApproachToString(TvdApproach::li) << ": li's approach for nonuniform cell sizes\n"
159 << tvdApproachToString(TvdApproach::hou) << ": hou's approach for nonuniform cell sizes");
160 }
161
166 {
167 switch (tvd)
168 {
169 case TvdApproach::uniform: return "Uniform";
170 case TvdApproach::li: return "Li";
171 case TvdApproach::hou: return "Hou";
172 default: return "Invalid"; // should never be reached
173 }
174 }
175
181 {
182 if (differencingScheme == "Vanleer") return DifferencingScheme::vanleer;
183 if (differencingScheme == "Vanalbada") return DifferencingScheme::vanalbada;
184 if (differencingScheme == "Minmod") return DifferencingScheme::minmod;
185 if (differencingScheme == "Superbee") return DifferencingScheme::superbee;
186 if (differencingScheme == "Umist") return DifferencingScheme::umist;
187 if (differencingScheme == "Mclimiter") return DifferencingScheme::mclimiter;
188 if (differencingScheme == "Wahyd") return DifferencingScheme::wahyd;
189 DUNE_THROW(ParameterException, "\nThis differencing scheme: \"" << differencingScheme << "\" is not implemented.\n"
190 << "The available differencing schemes are as follows: \n"
191 << differencingSchemeToString(DifferencingScheme::vanleer) << ": The Vanleer flux limiter\n"
192 << differencingSchemeToString(DifferencingScheme::vanalbada) << ": The VanAlbada flux limiter\n"
193 << differencingSchemeToString(DifferencingScheme::minmod) << ": The Min-Mod flux limiter\n"
194 << differencingSchemeToString(DifferencingScheme::superbee) << ": The SuperBEE flux limiter\n"
195 << differencingSchemeToString(DifferencingScheme::umist) << ": The UMist flux limiter\n"
196 << differencingSchemeToString(DifferencingScheme::mclimiter) << ": The McLimiter flux limiter\n"
197 << differencingSchemeToString(DifferencingScheme::wahyd) << ": The Wahyd flux limiter");
198 }
199
204 {
205 switch (differencingScheme)
206 {
207 case DifferencingScheme::vanleer: return "Vanleer";
208 case DifferencingScheme::vanalbada: return "Vanalbada";
209 case DifferencingScheme::minmod: return "Minmod";
210 case DifferencingScheme::superbee: return "Superbee";
211 case DifferencingScheme::umist: return "Umist";
212 case DifferencingScheme::mclimiter: return "Mclimiter";
213 case DifferencingScheme::wahyd: return "Wahyd";
214 default: return "Invalid"; // should never be reached
215 }
216 }
217
221 Scalar upwind(const Scalar downstreamMomentum,
222 const Scalar upstreamMomentum) const
223 {
224 return (upwindWeight_ * upstreamMomentum + (1.0 - upwindWeight_) * downstreamMomentum);
225 }
226
231 Scalar tvd(const std::array<Scalar,3>& momenta,
232 const std::array<Scalar,3>& distances,
233 const bool selfIsUpstream,
234 const TvdApproach tvdApproach) const
235 {
236 Scalar momentum = 0.0;
237 switch(tvdApproach)
238 {
240 {
241 momentum += tvdUniform(momenta, distances, selfIsUpstream);
242 break;
243 }
244 case TvdApproach::li :
245 {
246 momentum += tvdLi(momenta, distances, selfIsUpstream);
247 break;
248 }
249 case TvdApproach::hou :
250 {
251 momentum += tvdHou(momenta, distances, selfIsUpstream);
252 break;
253 }
254 default:
255 {
256 DUNE_THROW(ParameterException, "\nThis Tvd Approach is not implemented.\n");
257 break;
258 }
259 }
260 return momentum;
261 }
262
268 Scalar tvdUniform(const std::array<Scalar,3>& momenta,
269 const std::array<Scalar,3>& distances,
270 const bool selfIsUpstream) const
271 {
272 using std::isfinite;
273 const Scalar downstreamMomentum = momenta[0];
274 const Scalar upstreamMomentum = momenta[1];
275 const Scalar upUpstreamMomentum = momenta[2];
276 const Scalar ratio = (upstreamMomentum - upUpstreamMomentum) / (downstreamMomentum - upstreamMomentum);
277
278 // If the velocity field is uniform (like at the first newton step) we get a NaN
279 if(ratio > 0.0 && isfinite(ratio))
280 {
281 const Scalar secondOrderTerm = 0.5 * limiter_(ratio, 2.0) * (downstreamMomentum - upstreamMomentum);
282 return (upstreamMomentum + secondOrderTerm);
283 }
284 else
285 return upstreamMomentum;
286 }
287
295 Scalar tvdLi(const std::array<Scalar,3>& momenta,
296 const std::array<Scalar,3>& distances,
297 const bool selfIsUpstream) const
298 {
299 using std::isfinite;
300 const Scalar downstreamMomentum = momenta[0];
301 const Scalar upstreamMomentum = momenta[1];
302 const Scalar upUpstreamMomentum = momenta[2];
303 const Scalar upstreamToDownstreamDistance = distances[0];
304 const Scalar upUpstreamToUpstreamDistance = distances[1];
305
306 // selfIsUpstream is required to get the correct sign because upUpstreamToUpstreamDistance is always positive
307 const Scalar upUpstreamGradient = (upstreamMomentum - upUpstreamMomentum) / upUpstreamToUpstreamDistance * selfIsUpstream;
308
309 // Distance between the upUpstream node and the position where it should be if the grid were uniform.
310 const Scalar correctionDistance = upUpstreamToUpstreamDistance - upstreamToDownstreamDistance;
311 const Scalar reconstrutedUpUpstreamVelocity = upUpstreamMomentum + upUpstreamGradient * correctionDistance;
312 const Scalar ratio = (upstreamMomentum - reconstrutedUpUpstreamVelocity) / (downstreamMomentum - upstreamMomentum);
313
314 // If the velocity field is uniform (like at the first newton step) we get a NaN
315 if(ratio > 0.0 && isfinite(ratio))
316 {
317 const Scalar secondOrderTerm = 0.5 * limiter_(ratio, 2.0) * (downstreamMomentum - upstreamMomentum);
318 return (upstreamMomentum + secondOrderTerm);
319 }
320 else
321 return upstreamMomentum;
322 }
323
330 Scalar tvdHou(const std::array<Scalar,3>& momenta,
331 const std::array<Scalar,3>& distances,
332 const bool selfIsUpstream) const
333 {
334 using std::isfinite;
335 const Scalar downstreamMomentum = momenta[0];
336 const Scalar upstreamMomentum = momenta[1];
337 const Scalar upUpstreamMomentum = momenta[2];
338 const Scalar upstreamToDownstreamDistance = distances[0];
339 const Scalar upUpstreamToUpstreamDistance = distances[1];
340 const Scalar downstreamStaggeredCellSize = distances[2];
341 const Scalar ratio = (upstreamMomentum - upUpstreamMomentum) / (downstreamMomentum - upstreamMomentum)
342 * upstreamToDownstreamDistance / upUpstreamToUpstreamDistance;
343
344 // If the velocity field is uniform (like at the first newton step) we get a NaN
345 if(ratio > 0.0 && isfinite(ratio))
346 {
347 const Scalar upstreamStaggeredCellSize = 0.5 * (upstreamToDownstreamDistance + upUpstreamToUpstreamDistance);
348 const Scalar R = (upstreamStaggeredCellSize + downstreamStaggeredCellSize) / upstreamStaggeredCellSize;
349 const Scalar secondOrderTerm = limiter_(ratio, R) / R * (downstreamMomentum - upstreamMomentum);
350 return (upstreamMomentum + secondOrderTerm);
351 }
352 else
353 return upstreamMomentum;
354 }
355
361 static Scalar vanleer(const Scalar r, const Scalar R)
362 {
363 return R * r / (R - 1.0 + r);
364 }
365
369 static Scalar vanalbada(const Scalar r, const Scalar R)
370 {
371 return r * (r + 1.0) / (1.0 + r * r);
372 }
373
377 static Scalar minmod(const Scalar r, const Scalar R)
378 {
379 using std::min;
380 return min(r, 1.0);
381 }
382
388 static Scalar superbee(const Scalar r, const Scalar R)
389 {
390 using std::min;
391 using std::max;
392 return max(min(R * r, 1.0), min(r, R));
393 }
394
398 static Scalar umist(const Scalar r, const Scalar R)
399 {
400 using std::min;
401 return min({R * r, (r * (5.0 - R) + R - 1.0) / 4.0, (r * (R - 1.0) + 5.0 - R) / 4.0, R});
402 }
403
404 /*
405 * \brief Monotonized-Central limiter [Van Leer 1977]
406 */
407 static Scalar mclimiter(const Scalar r, const Scalar R)
408 {
409 using std::min;
410 return min({R * r, (r + 1.0) / 2.0, R});
411 }
412
416 static Scalar wahyd(const Scalar r, const Scalar R)
417 {
418 using std::min;
419 return r > 1 ? min((r + R * r * r) / (R + r * r), R)
420 : vanleer(r, R);
421 }
422
425 {
426 return tvdApproach_;
427 }
428
431 {
432 return differencingScheme_;
433 }
434
435private:
436 TvdApproach tvdApproach_;
437 DifferencingScheme differencingScheme_;
438 Scalar upwindWeight_;
439
440 std::function<Scalar(const Scalar, const Scalar)> limiter_;
441};
442
443} // end namespace Dumux
444
445#endif
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:48
This file contains different higher order methods for approximating the velocity.
Definition: staggeredupwindmethods.hh:50
TvdApproach tvdApproachFromString(const std::string &tvd)
Convenience function to convert user input given as std::string to the corresponding enum class used ...
Definition: staggeredupwindmethods.hh:150
static Scalar vanalbada(const Scalar r, const Scalar R)
Van Albada flux limiter function [Van Albada et al. 1982].
Definition: staggeredupwindmethods.hh:369
StaggeredUpwindMethods(const std::string &paramGroup="")
Definition: staggeredupwindmethods.hh:52
Scalar upwind(const Scalar downstreamMomentum, const Scalar upstreamMomentum) const
Upwind Method.
Definition: staggeredupwindmethods.hh:221
const DifferencingScheme & differencingScheme() const
Returns the differencing scheme.
Definition: staggeredupwindmethods.hh:430
std::string tvdApproachToString(TvdApproach tvd)
return the name of the TVD approach
Definition: staggeredupwindmethods.hh:165
static Scalar umist(const Scalar r, const Scalar R)
UMIST flux limiter function [Lien and Leschziner 1993].
Definition: staggeredupwindmethods.hh:398
Scalar tvdUniform(const std::array< Scalar, 3 > &momenta, const std::array< Scalar, 3 > &distances, const bool selfIsUpstream) const
Tvd Scheme: Total Variation Diminishing.
Definition: staggeredupwindmethods.hh:268
Scalar tvdLi(const std::array< Scalar, 3 > &momenta, const std::array< Scalar, 3 > &distances, const bool selfIsUpstream) const
Tvd Scheme: Total Variation Diminishing.
Definition: staggeredupwindmethods.hh:295
static Scalar mclimiter(const Scalar r, const Scalar R)
Definition: staggeredupwindmethods.hh:407
DifferencingScheme differencingSchemeFromString(const std::string &differencingScheme)
Convenience function to convert user input given as std::string to the corresponding enum class used ...
Definition: staggeredupwindmethods.hh:180
Scalar tvd(const std::array< Scalar, 3 > &momenta, const std::array< Scalar, 3 > &distances, const bool selfIsUpstream, const TvdApproach tvdApproach) const
Tvd Scheme: Total Variation Diminishing.
Definition: staggeredupwindmethods.hh:231
Scalar tvdHou(const std::array< Scalar, 3 > &momenta, const std::array< Scalar, 3 > &distances, const bool selfIsUpstream) const
Tvd Scheme: Total Variation Diminishing.
Definition: staggeredupwindmethods.hh:330
static Scalar vanleer(const Scalar r, const Scalar R)
Van Leer flux limiter function [Van Leer 1974].
Definition: staggeredupwindmethods.hh:361
const TvdApproach & tvdApproach() const
Returns the Tvd approach.
Definition: staggeredupwindmethods.hh:424
std::string differencingSchemeToString(DifferencingScheme differencingScheme)
return the name of the Discretization Method
Definition: staggeredupwindmethods.hh:203
static Scalar superbee(const Scalar r, const Scalar R)
SUPERBEE flux limiter function [Roe 1985].
Definition: staggeredupwindmethods.hh:388
static Scalar minmod(const Scalar r, const Scalar R)
MinMod flux limiter function [Roe 1985].
Definition: staggeredupwindmethods.hh:377
static Scalar wahyd(const Scalar r, const Scalar R)
WAHYD Scheme [Hou, Simons, Hinkelmann 2007];.
Definition: staggeredupwindmethods.hh:416
Some exceptions thrown in DuMux
TvdApproach
Available Tvd approaches.
Definition: staggeredupwindmethods.hh:31
DifferencingScheme
Available differencing schemes.
Definition: staggeredupwindmethods.hh:40
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.