25#ifndef DUMUX_UPWINDING_METHODS_HH
26#define DUMUX_UPWINDING_METHODS_HH
60template<
class Scalar,
int upwindSchemeOrder>
66 upwindWeight_ = getParamFromGroup<Scalar>(paramGroup,
"Flux.UpwindWeight");
68 if (upwindSchemeOrder > 1)
71 tvdApproach_ =
tvdApproachFromString(getParamFromGroup<std::string>(paramGroup,
"Flux.TvdApproach",
"Uniform"));
75 switch (differencingScheme_)
85 DUNE_THROW(
ParameterException,
"\nDifferencing scheme (Van Albada) is not implemented for the Tvd approach (Hou).");
102 limiter_ = this->
umist;
112 limiter_ = this->
wahyd;
118 " is not implemented.\n");
125 std::cout <<
"No TvdApproach specified. Defaulting to the Uniform method." <<
"\n";
126 std::cout <<
"Other available TVD approaches for uniform (and nonuniform) grids are as follows: \n"
130 std::cout <<
"Each approach can be specified as written above in the Flux group under the title TvdApproach in your input file. \n";
135 std::cout <<
"No DifferencingScheme specified. Defaulting to the Minmod scheme." <<
"\n";
136 std::cout <<
"Other available Differencing Schemes are as follows: \n"
144 std::cout <<
"Each scheme can be specified as written above in the Flux group under the variable DifferencingScheme in your input file. \n";
168 <<
"The available TVD approaches for uniform (and nonuniform) grids are as follows: \n"
184 default:
return "Invalid";
202 <<
"The available differencing schemes are as follows: \n"
226 default:
return "Invalid";
233 Scalar
upwind(
const Scalar downstreamMomentum,
234 const Scalar upstreamMomentum)
const
236 return (upwindWeight_ * upstreamMomentum + (1.0 - upwindWeight_) * downstreamMomentum);
243 Scalar
tvd(
const std::array<Scalar,3>& momenta,
244 const std::array<Scalar,3>& distances,
245 const bool selfIsUpstream,
248 Scalar momentum = 0.0;
253 momentum +=
tvdUniform(momenta, distances, selfIsUpstream);
258 momentum +=
tvdLi(momenta, distances, selfIsUpstream);
263 momentum +=
tvdHou(momenta, distances, selfIsUpstream);
281 const std::array<Scalar,3>& distances,
282 const bool selfIsUpstream)
const
285 const Scalar downstreamMomentum = momenta[0];
286 const Scalar upstreamMomentum = momenta[1];
287 const Scalar upUpstreamMomentum = momenta[2];
288 const Scalar ratio = (upstreamMomentum - upUpstreamMomentum) / (downstreamMomentum - upstreamMomentum);
291 if(ratio > 0.0 && isfinite(ratio))
293 const Scalar secondOrderTerm = 0.5 * limiter_(ratio, 2.0) * (downstreamMomentum - upstreamMomentum);
294 return (upstreamMomentum + secondOrderTerm);
297 return upstreamMomentum;
307 Scalar
tvdLi(
const std::array<Scalar,3>& momenta,
308 const std::array<Scalar,3>& distances,
309 const bool selfIsUpstream)
const
312 const Scalar downstreamMomentum = momenta[0];
313 const Scalar upstreamMomentum = momenta[1];
314 const Scalar upUpstreamMomentum = momenta[2];
315 const Scalar upstreamToDownstreamDistance = distances[0];
316 const Scalar upUpstreamToUpstreamDistance = distances[1];
319 const Scalar upUpstreamGradient = (upstreamMomentum - upUpstreamMomentum) / upUpstreamToUpstreamDistance * selfIsUpstream;
322 const Scalar correctionDistance = upUpstreamToUpstreamDistance - upstreamToDownstreamDistance;
323 const Scalar reconstrutedUpUpstreamVelocity = upUpstreamMomentum + upUpstreamGradient * correctionDistance;
324 const Scalar ratio = (upstreamMomentum - reconstrutedUpUpstreamVelocity) / (downstreamMomentum - upstreamMomentum);
327 if(ratio > 0.0 && isfinite(ratio))
329 const Scalar secondOrderTerm = 0.5 * limiter_(ratio, 2.0) * (downstreamMomentum - upstreamMomentum);
330 return (upstreamMomentum + secondOrderTerm);
333 return upstreamMomentum;
342 Scalar
tvdHou(
const std::array<Scalar,3>& momenta,
343 const std::array<Scalar,3>& distances,
344 const bool selfIsUpstream)
const
347 const Scalar downstreamMomentum = momenta[0];
348 const Scalar upstreamMomentum = momenta[1];
349 const Scalar upUpstreamMomentum = momenta[2];
350 const Scalar upstreamToDownstreamDistance = distances[0];
351 const Scalar upUpstreamToUpstreamDistance = distances[1];
352 const Scalar downstreamStaggeredCellSize = distances[2];
353 const Scalar ratio = (upstreamMomentum - upUpstreamMomentum) / (downstreamMomentum - upstreamMomentum)
354 * upstreamToDownstreamDistance / upUpstreamToUpstreamDistance;
357 if(ratio > 0.0 && isfinite(ratio))
359 const Scalar upstreamStaggeredCellSize = 0.5 * (upstreamToDownstreamDistance + upUpstreamToUpstreamDistance);
360 const Scalar R = (upstreamStaggeredCellSize + downstreamStaggeredCellSize) / upstreamStaggeredCellSize;
361 const Scalar secondOrderTerm = limiter_(ratio, R) / R * (downstreamMomentum - upstreamMomentum);
362 return (upstreamMomentum + secondOrderTerm);
365 return upstreamMomentum;
373 static Scalar
vanleer(
const Scalar r,
const Scalar R)
375 return R * r / (R - 1.0 + r);
383 return r * (r + 1.0) / (1.0 + r * r);
389 static Scalar
minmod(
const Scalar r,
const Scalar R)
400 static Scalar
superbee(
const Scalar r,
const Scalar R)
404 return max(min(R * r, 1.0), min(r, R));
410 static Scalar
umist(
const Scalar r,
const Scalar R)
413 return min({R * r, (r * (5.0 - R) + R - 1.0) / 4.0, (r * (R - 1.0) + 5.0 - R) / 4.0, R});
422 return min({R * r, (r + 1.0) / 2.0, R});
428 static Scalar
wahyd(
const Scalar r,
const Scalar R)
431 return r > 1 ? min((r + R * r * r) / (R + r * r), R)
444 return differencingScheme_;
450 Scalar upwindWeight_;
452 std::function<Scalar(
const Scalar,
const Scalar)> limiter_;
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
TvdApproach
Available Tvd approaches.
Definition: staggeredupwindmethods.hh:43
DifferencingScheme
Available differencing schemes.
Definition: staggeredupwindmethods.hh:52
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:390
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60
This file contains different higher order methods for approximating the velocity.
Definition: staggeredupwindmethods.hh:62
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:162
static Scalar vanalbada(const Scalar r, const Scalar R)
Van Albada flux limiter function [Van Albada et al. 1982].
Definition: staggeredupwindmethods.hh:381
StaggeredUpwindMethods(const std::string ¶mGroup="")
Definition: staggeredupwindmethods.hh:64
Scalar upwind(const Scalar downstreamMomentum, const Scalar upstreamMomentum) const
Upwind Method.
Definition: staggeredupwindmethods.hh:233
const DifferencingScheme & differencingScheme() const
Returns the differencing scheme.
Definition: staggeredupwindmethods.hh:442
std::string tvdApproachToString(TvdApproach tvd)
return the name of the TVD approach
Definition: staggeredupwindmethods.hh:177
static Scalar umist(const Scalar r, const Scalar R)
UMIST flux limiter function [Lien and Leschziner 1993].
Definition: staggeredupwindmethods.hh:410
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:280
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:307
static Scalar mclimiter(const Scalar r, const Scalar R)
Definition: staggeredupwindmethods.hh:419
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:192
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:243
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:342
static Scalar vanleer(const Scalar r, const Scalar R)
Van Leer flux limiter function [Van Leer 1974].
Definition: staggeredupwindmethods.hh:373
const TvdApproach & tvdApproach() const
Returns the Tvd approach.
Definition: staggeredupwindmethods.hh:436
std::string differencingSchemeToString(DifferencingScheme differencingScheme)
return the name of the Discretization Method
Definition: staggeredupwindmethods.hh:215
static Scalar superbee(const Scalar r, const Scalar R)
SUPERBEE flux limiter function [Roe 1985].
Definition: staggeredupwindmethods.hh:400
static Scalar minmod(const Scalar r, const Scalar R)
MinMod flux limiter function [Roe 1985].
Definition: staggeredupwindmethods.hh:389
static Scalar wahyd(const Scalar r, const Scalar R)
WAHYD Scheme [Hou, Simons, Hinkelmann 2007];.
Definition: staggeredupwindmethods.hh:428