3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Files | Functions
Flux related to the shallow water model

Flux related to the shallow water model. More...

Description

Flux related to the shallow water model.

Files

file  exactriemann.hh
 Function to compute the Riemann flux at the interface.
 
file  fluxlimiterlet.hh
 Function to limit the fluxes.
 
file  riemannproblem.hh
 This file includes a function to construct the Riemann problem.
 

Functions

template<class Scalar >
RiemannSolution< Scalar > Dumux::ShallowWater::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. More...
 
template<class Scalar >
static Scalar Dumux::ShallowWater::fluxLimiterLET (const Scalar valueLeft, const Scalar valueRight, const Scalar upperH, const Scalar lowerH)
 Flux limiter function to scale fluxes for small water depths. More...
 
template<class Scalar , class GlobalPosition >
std::array< Scalar, 3 > Dumux::ShallowWater::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. More...
 

Function Documentation

◆ exactRiemann()

template<class Scalar >
RiemannSolution< Scalar > Dumux::ShallowWater::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.

The flux of the 2D shallow water equations must be rotated to a 1D problem before the Riemann solver can be applied. The computed water flux is given in m^2/s, the momentum fluxes are given in m^3/s^2.

This Riemann solver is described in the book "Shock-capturing methods for free-surface shallow flows" from Toro, 2001. We keep the notation for the variables after Toro.

Parameters
dlwater depth on the left side
drwater depth on the right side
ulveloctiyX on the left side
urvelocityX on the right side
vlvelocityY on the left side
vrvelocityY on the right side
gravgravity constant
ssample point (default = 0 since x = 0 for flux computation)

◆ fluxLimiterLET()

template<class Scalar >
static Scalar Dumux::ShallowWater::fluxLimiterLET ( const Scalar  valueLeft,
const Scalar  valueRight,
const Scalar  upperH,
const Scalar  lowerH 
)
static

Flux limiter function to scale fluxes for small water depths.

This function acts like a kind of mobility, it limits the water flux for small water depths. The mobility depends on the left and right side state of a variable. The LET-Type function is described at https://en.wikipedia.org/wiki/Relative_permeability The LET-Parameters are fixed. The current parameters are set to a values to get a nice curve. They have no physical meaning.

Template Parameters
Scalarthe scalar type for scalar physical quantities
Parameters
valueLeftThe value on the left side
valueRightThe value on the right side
upperHWhere to start the limit the function (mobility < 1)
lowerHWhere the limit should reach zero (mobility < 0)

◆ riemannProblem()

template<class Scalar , class GlobalPosition >
std::array< Scalar, 3 > Dumux::ShallowWater::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.

Riemann problem applies the hydrostatic reconstruction, uses the Riemann invariants to transform the two-dimensional problem to a one-dimensional problem, solves this new problem, and rotates the problem back. Further it applies a flux limiter for the water flux to handle drying elements. The correction of the bed slope source term leads to a non-symmetric flux term at the interface for the momentum equations. Since DuMuX computes the fluxes twice from each side this does not matter.

So far this implements the exact Riemann solver (with reconstruction after Audusse).

The computed water flux (localFlux[0]) is given in m^2/s, the momentum fluxes (localFlux[1], localFlux[2]) are given in m^3/s^2. Later this flux will be multiplied by the scvf.area() (given in m for a 2D problem) to get the flux over a face.

Parameters
waterDepthLeftwater depth on the left side
waterDepthRightwater depth on the right side
velocityXLeftveloctiyX on the left side
velocityXRightvelocityX on the right side
velocityYLeftvelocityY on the left side
velocityYRightvelocityY on the right side
bedSurfaceLeftsurface of the bed on the left side
bedSurfaceRightsurface of the bed on the right side
gravitygravity constant
nxythe normal vector