Linear Shallow Water Equations#
Definition#
The SELF_LinearShallowWater2D_t
module defines the LinearShallowWater2D_t
class. In SELF, models are posed in the form of a generic conservation law
where \(\vec{s}\) is a vector of solution variables, \(\overleftrightarrow{f}\) is the conservative flux, and \(\vec{q}\) are non-conservative source terms.
For the linear shallow water equations in 2-D
where \(u\) and \(v\) are the x and y components of the barotropic velocity (\(\vec{u} = u \hat{x} + v \hat{y}\)) and \(\eta\) is the deviation of the fluid free surface relative to the resting fluid.
where \(g\) is acceleration due to gravity and \(H\) is uniform resting fluid depth. The source term is set to zero.
To track stability of the Euler equation, the total entropy function is
Implementation#
The 2D Linear Shallow Water model is implemented as a type extension of the DGModel2d
class. The LinearShallowWater2D_t
class adds parameters for acceleration due to gravity and the uniform resting fluid depth. It also overrides SetMetaData
, entropy_func
, flux2d
, and riemannflux2d
type-bound procedures.
Riemann Solver#
The LinearShallowWater2D
class is defined using the advective form.
The Riemann solver for the hyperbolic part of the shallow water equations is the local Lax-Friedrichs upwind Riemann solver
where \(c = \sqrt{gH}\), and
Together, this becomes
The details for this implementation can be found in self_LinearShallowWater2D_t.f90.
Boundary Conditions#
When initializing the mesh for your 2D Linear Shallow Water Equations solver, you can change the boundary conditions to
SELF_BC_Radiation
to set the external state on model boundaries to 0 in the Riemann solverSELF_BC_NoNormalFlow
to set the external normal velocity to the negative of the interior normal velocity and prolong the density, pressure, and tangential velocity (free slip). This effectively creates a reflecting boundary condition.SELF_BC_Prescribed
to set a prescribed external state.
As an example, when using the built-in structured mesh generator, you can do the following
type(Mesh2D),target :: mesh
integer :: bcids(1:4)
bcids(1:4) = (/&
SELF_NONORMALFLOW,& ! South boundary condition
SELF_RADIATION,& ! East boundary condition
SELF_PRESCRIBED,& ! North boundary condition
SELF_RADIATION & ! West boundary condition
/)
call mesh%StructuredMesh(nxPerTile=5,nyPerTile=5,&
nTileX=2,nTileY=2,&
dx=0.1_prec,dy=0.1_prec,bcids)
Note
See the Structured Mesh documentation for details on using the structuredmesh
procedure
Note
To set a prescribed state as a function of position and time, you can create a type-extension of the LinearShallowWater2D
class and override the hbc2d_Prescribed
Example usage#
For examples, see any of the following
examples/LinearShallowWater2D.f90
- Implements the 2D shallow water equations with no normal flow boundary conditions