3.5-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
70template<class Scalar, class GlobalPosition>
71std::array<Scalar,3> riemannProblem(const Scalar waterDepthLeft,
72 const Scalar waterDepthRight,
73 Scalar velocityXLeft,
74 Scalar velocityXRight,
75 Scalar velocityYLeft,
76 Scalar velocityYRight,
77 const Scalar bedSurfaceLeft,
78 const Scalar bedSurfaceRight,
79 const Scalar gravity,
80 const GlobalPosition& nxy)
81{
82 using std::max;
83
84 // hydrostatic reconstrucion after Audusse
85 const Scalar dzl = max(0.0, bedSurfaceRight - bedSurfaceLeft);
86 const Scalar waterDepthLeftReconstructed = max(0.0, waterDepthLeft - dzl);
87 const Scalar dzr = max(0.0, bedSurfaceLeft - bedSurfaceRight);
88 const Scalar waterDepthRightReconstructed = max(0.0, waterDepthRight - dzr);
89
90 // make rotation of the flux we compute an 1d flux
91 Scalar tempFlux = velocityXLeft;
92 velocityXLeft = nxy[0] * tempFlux + nxy[1] * velocityYLeft;
93 velocityYLeft = -nxy[1] * tempFlux + nxy[0] * velocityYLeft;
94
95 tempFlux = velocityXRight;
96 velocityXRight = nxy[0] * tempFlux + nxy[1] * velocityYRight;
97 velocityYRight = -nxy[1] * tempFlux + nxy[0] * velocityYRight;
98
99 auto riemannResult = ShallowWater::exactRiemann(waterDepthLeftReconstructed,
100 waterDepthRightReconstructed,
101 velocityXLeft,
102 velocityXRight,
103 velocityYLeft,
104 velocityYRight,
105 gravity);
106
107 //redo rotation
108 tempFlux = riemannResult.flux[1];
109 riemannResult.flux[1] = nxy[0] * tempFlux - nxy[1] * riemannResult.flux[2];
110 riemannResult.flux[2] = nxy[1] * tempFlux + nxy[0] * riemannResult.flux[2];
111
112 // Add reconstruction flux from Audusse reconstruction
113 const Scalar hgzl = 0.5 * (waterDepthLeftReconstructed + waterDepthLeft) * (waterDepthLeftReconstructed - waterDepthLeft);
114 const Scalar hdxzl = gravity * nxy[0] * hgzl;
115 const Scalar hdyzl = gravity * nxy[1] * hgzl;
116
117 /*Right side is computed from the other side otherwise the
118 following "non-symetric" fluxes are needed:
119
120 Scalar hgzr = 0.5 * (waterDepthRightReconstructed + waterDepthRight) * (waterDepthRightReconstructed - waterDepthRight);
121 Scalar hdxzr = gravity * nxy[0] * hgzr;
122 Scalar hdyzrhdyzr = gravity * nxy[1] * hgzr;
123 */
124
125 // compute the mobility of the flux with the fluxlimiter
126 static const Scalar upperWaterDepthFluxLimiting = getParam<Scalar>("FluxLimiterLET.UpperWaterDepth", 1e-3);
127 static const Scalar lowerWaterDepthFluxLimiting = getParam<Scalar>("FluxLimiterLET.LowerWaterDepth", 1e-5);
128 static const bool upwindWaterDepthFluxLimiting = getParam<bool>("FluxLimiterLET.UpwindFluxLimiting", false);
129
130 Scalar limitingDepth = (waterDepthLeftReconstructed + waterDepthRightReconstructed) * 0.5;
131
132 //Using the upwind water depth from the flux direction can improve stability.
133 if (upwindWaterDepthFluxLimiting)
134 {
135 if (riemannResult.flux[0] < 0)
136 {
137 limitingDepth = waterDepthRightReconstructed;
138 }else
139 {
140 limitingDepth = waterDepthLeftReconstructed;
141 }
142 }
143
144 const Scalar mobility = ShallowWater::fluxLimiterLET(limitingDepth,
145 limitingDepth,
146 upperWaterDepthFluxLimiting,
147 lowerWaterDepthFluxLimiting);
148 std::array<Scalar, 3> localFlux;
149 localFlux[0] = riemannResult.flux[0] * mobility;
150 localFlux[1] = (riemannResult.flux[1] - hdxzl);
151 localFlux[2] = (riemannResult.flux[2] - hdyzl);
152
153 return localFlux;
154}
155
156} // end namespace ShallowWater
157} // end namespace Dumux
158
159#endif
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:71
Definition: adapt.hh:29
std::string mobility(int phaseIdx) noexcept
I/O name of mobility for multiphase systems.
Definition: name.hh:101