25#ifndef DUMUX_FLUX_SHALLOW_WATER_RIEMANN_PROBLEM_HH
26#define DUMUX_FLUX_SHALLOW_WATER_RIEMANN_PROBLEM_HH
33namespace ShallowWater {
71template<
class Scalar,
class GlobalPosition>
73 const Scalar waterDepthRight,
75 Scalar velocityXRight,
77 Scalar velocityYRight,
78 const Scalar bedSurfaceLeft,
79 const Scalar bedSurfaceRight,
81 const GlobalPosition& nxy)
86 const Scalar dzl = max(0.0, bedSurfaceRight - bedSurfaceLeft);
87 const Scalar waterDepthLeftReconstructed = max(0.0, waterDepthLeft - dzl);
88 const Scalar dzr = max(0.0, bedSurfaceLeft - bedSurfaceRight);
89 const Scalar waterDepthRightReconstructed = max(0.0, waterDepthRight - dzr);
92 Scalar tempFlux = velocityXLeft;
93 velocityXLeft = nxy[0] * tempFlux + nxy[1] * velocityYLeft;
94 velocityYLeft = -nxy[1] * tempFlux + nxy[0] * velocityYLeft;
96 tempFlux = velocityXRight;
97 velocityXRight = nxy[0] * tempFlux + nxy[1] * velocityYRight;
98 velocityYRight = -nxy[1] * tempFlux + nxy[0] * velocityYRight;
101 waterDepthRightReconstructed,
109 tempFlux = riemannResult.flux[1];
110 riemannResult.flux[1] = nxy[0] * tempFlux - nxy[1] * riemannResult.flux[2];
111 riemannResult.flux[2] = nxy[1] * tempFlux + nxy[0] * riemannResult.flux[2];
114 const Scalar hgzl = 0.5 * (waterDepthLeftReconstructed + waterDepthLeft) * (waterDepthLeftReconstructed - waterDepthLeft);
115 const Scalar hdxzl = gravity * nxy[0] * hgzl;
116 const Scalar hdyzl = gravity * nxy[1] * hgzl;
127 static const Scalar upperWaterDepthFluxLimiting = getParam<Scalar>(
"FluxLimiterLET.UpperWaterDepth", 1e-3);
128 static const Scalar lowerWaterDepthFluxLimiting = getParam<Scalar>(
"FluxLimiterLET.LowerWaterDepth", 1e-5);
129 static const bool upwindWaterDepthFluxLimiting = getParam<bool>(
"FluxLimiterLET.UpwindFluxLimiting",
false);
131 Scalar limitingDepth = (waterDepthLeftReconstructed + waterDepthRightReconstructed) * 0.5;
134 if (upwindWaterDepthFluxLimiting)
136 if (riemannResult.flux[0] < 0)
138 limitingDepth = waterDepthRightReconstructed;
141 limitingDepth = waterDepthLeftReconstructed;
147 upperWaterDepthFluxLimiting,
148 lowerWaterDepthFluxLimiting);
149 std::array<Scalar, 3> localFlux;
150 localFlux[0] = riemannResult.flux[0] *
mobility;
151 localFlux[1] = (riemannResult.flux[1] - hdxzl);
152 localFlux[2] = (riemannResult.flux[2] - hdyzl);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Function to limit the fluxes.
Function to compute the Riemann flux at the interface.
RiemannSolution< Scalar > exactRiemann(const Scalar dl, const Scalar dr, const Scalar ul, const Scalar ur, const Scalar vl, const Scalar vr, const Scalar grav, const Scalar s=0.0)
Exact Riemann solver for the shallow water equations.
Definition: exactriemann.hh:66
static Scalar fluxLimiterLET(const Scalar valueLeft, const Scalar valueRight, const Scalar upperH, const Scalar lowerH)
Flux limiter function to scale fluxes for small water depths.
Definition: fluxlimiterlet.hh:52
std::array< Scalar, 3 > riemannProblem(const Scalar waterDepthLeft, const Scalar waterDepthRight, Scalar velocityXLeft, Scalar velocityXRight, Scalar velocityYLeft, Scalar velocityYRight, const Scalar bedSurfaceLeft, const Scalar bedSurfaceRight, const Scalar gravity, const GlobalPosition &nxy)
Construct a Riemann problem and solve it.
Definition: riemannproblem.hh:72
std::string mobility(int phaseIdx) noexcept
I/O name of mobility for multiphase systems.
Definition: name.hh:101