3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
riemannproblem.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 2 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_FLUX_SHALLOW_WATER_RIEMANN_PROBLEM_HH
26#define DUMUX_FLUX_SHALLOW_WATER_RIEMANN_PROBLEM_HH
27
31
32namespace Dumux {
33namespace ShallowWater {
34
66template<class Scalar, class GlobalPosition>
67std::array<Scalar,3> riemannProblem(const Scalar waterDepthLeft,
68 const Scalar waterDepthRight,
69 Scalar velocityXLeft,
70 Scalar velocityXRight,
71 Scalar velocityYLeft,
72 Scalar velocityYRight,
73 const Scalar bedSurfaceLeft,
74 const Scalar bedSurfaceRight,
75 const Scalar gravity,
76 const GlobalPosition& nxy)
77{
78 using std::max;
79
80 // hydrostatic reconstrucion after Audusse
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);
85
86 // make rotation of the flux we compute an 1d flux
87 Scalar tempFlux = velocityXLeft;
88 velocityXLeft = nxy[0] * tempFlux + nxy[1] * velocityYLeft;
89 velocityYLeft = -nxy[1] * tempFlux + nxy[0] * velocityYLeft;
90
91 tempFlux = velocityXRight;
92 velocityXRight = nxy[0] * tempFlux + nxy[1] * velocityYRight;
93 velocityYRight = -nxy[1] * tempFlux + nxy[0] * velocityYRight;
94
95 auto riemannResult = ShallowWater::exactRiemann(waterDepthLeftReconstructed,
96 waterDepthRightReconstructed,
97 velocityXLeft,
98 velocityXRight,
99 velocityYLeft,
100 velocityYRight,
101 gravity);
102
103 //redo rotation
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];
107
108 // Add reconstruction flux from Audusse reconstruction
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;
112
113 /*Right side is computed from the other side otherwise the
114 following "non-symetric" fluxes are needed:
115
116 Scalar hgzr = 0.5 * (waterDepthRightReconstructed + waterDepthRight) * (waterDepthRightReconstructed - waterDepthRight);
117 Scalar hdxzr = gravity * nxy[0] * hgzr;
118 Scalar hdyzrhdyzr = gravity * nxy[1] * hgzr;
119 */
120
121 // compute the mobility of the flux with the fluxlimiter
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);
125
126 Scalar limitingDepth = (waterDepthLeftReconstructed + waterDepthRightReconstructed) * 0.5;
127
128 //Using the upwind water depth from the flux direction can improve stability.
129 if (upwindWaterDepthFluxLimiting)
130 {
131 if (riemannResult.flux[0] < 0)
132 {
133 limitingDepth = waterDepthRightReconstructed;
134 }else
135 {
136 limitingDepth = waterDepthLeftReconstructed;
137 }
138 }
139
140 const Scalar mobility = ShallowWater::fluxLimiterLET(limitingDepth,
141 limitingDepth,
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);
148
149 return localFlux;
150}
151
152} // end namespace ShallowWater
153} // end namespace Dumux
154
155#endif
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
Definition: adapt.hh:29
std::string mobility(int phaseIdx) noexcept
I/O name of mobility for multiphase systems.
Definition: name.hh:101