25#ifndef DUMUX_FLUX_SHALLOW_WATER_RIEMANN_PROBLEM_HH
26#define DUMUX_FLUX_SHALLOW_WATER_RIEMANN_PROBLEM_HH
33namespace ShallowWater {
66template<
class Scalar,
class GlobalPosition>
68 const Scalar waterDepthRight,
70 Scalar velocityXRight,
72 Scalar velocityYRight,
73 const Scalar bedSurfaceLeft,
74 const Scalar bedSurfaceRight,
76 const GlobalPosition& nxy)
81 const Scalar dzl = max(0.0, bedSurfaceRight - bedSurfaceLeft);
82 const Scalar waterDepthLeftReconstructed = max(0.0, waterDepthLeft - dzl);
83 const Scalar dzr = max(0.0, bedSurfaceLeft - bedSurfaceRight);
84 const Scalar waterDepthRightReconstructed = max(0.0, waterDepthRight - dzr);
87 Scalar tempFlux = velocityXLeft;
88 velocityXLeft = nxy[0] * tempFlux + nxy[1] * velocityYLeft;
89 velocityYLeft = -nxy[1] * tempFlux + nxy[0] * velocityYLeft;
91 tempFlux = velocityXRight;
92 velocityXRight = nxy[0] * tempFlux + nxy[1] * velocityYRight;
93 velocityYRight = -nxy[1] * tempFlux + nxy[0] * velocityYRight;
96 waterDepthRightReconstructed,
104 tempFlux = riemannResult.flux[1];
105 riemannResult.flux[1] = nxy[0] * tempFlux - nxy[1] * riemannResult.flux[2];
106 riemannResult.flux[2] = nxy[1] * tempFlux + nxy[0] * riemannResult.flux[2];
109 const Scalar hgzl = 0.5 * (waterDepthLeftReconstructed + waterDepthLeft) * (waterDepthLeftReconstructed - waterDepthLeft);
110 const Scalar hdxzl = gravity * nxy[0] * hgzl;
111 const Scalar hdyzl = gravity * nxy[1] * hgzl;
122 static const Scalar upperWaterDepthFluxLimiting = getParam<Scalar>(
"FluxLimiterLET.UpperWaterDepth", 1e-3);
123 static const Scalar lowerWaterDepthFluxLimiting = getParam<Scalar>(
"FluxLimiterLET.LowerWaterDepth", 1e-5);
124 static const bool upwindWaterDepthFluxLimiting = getParam<bool>(
"FluxLimiterLET.UpwindFluxLimiting",
false);
126 Scalar limitingDepth = (waterDepthLeftReconstructed + waterDepthRightReconstructed) * 0.5;
129 if (upwindWaterDepthFluxLimiting)
131 if (riemannResult.flux[0] < 0)
133 limitingDepth = waterDepthRightReconstructed;
136 limitingDepth = waterDepthLeftReconstructed;
142 upperWaterDepthFluxLimiting,
143 lowerWaterDepthFluxLimiting);
144 std::array<Scalar, 3> localFlux;
145 localFlux[0] = riemannResult.flux[0] *
mobility;
146 localFlux[1] = (riemannResult.flux[1] - hdxzl);
147 localFlux[2] = (riemannResult.flux[2] - hdyzl);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Function to compute the Riemann flux at the interface.
Function to limit the fluxes.
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 Shallow water equations.
Definition: exactriemann.hh:61
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 Riemann Problem and solve it.
Definition: riemannproblem.hh:67
std::string mobility(int phaseIdx) noexcept
I/O name of mobility for multiphase systems.
Definition: name.hh:101