Shallow water Riemann solvers in Clawpack

A wide range of shallow water (SW) solvers are available in clawpack.riemann. Here’s a brief description of each. For each one, we have indicated (after “Fortran:”) the files you should compile to use it in the Fortran codes, and after “PyClaw” where you should import it from to use it in Python. If a pure-Python implementation is available, we also indicate that. Finally, we include links to examples that use each solver.

One dimension

For most 1D solvers, the vector q of conserved quantities is

\[\begin{split}q = \begin{bmatrix} h \\ hu \end{bmatrix},\end{split}\]

where \(h\) is depth and \(hu\) is momentum. Solvers with a tracer include that as a 3rd component. For solvers with bathymetry, the bathymetry is the first (and only) component of aux. All solvers require setting a constant parameter grav, which controls the force of gravity.

Two dimensions

For most 2D solvers, the vector q of conserved quantities is

\[\begin{split}q = \begin{bmatrix} h \\ hu \\ hv \end{bmatrix},\end{split}\]

where \(h\) is depth and \(hu, hv\) are the \(x\)- and \(y\)-components of momentum. For solvers with bathymetry, the bathymetry is the first (and only) component of aux. For the mapped solver, see the implementation for a description of aux. As in 1D, all solvers require setting a constant parameter grav, which controls the force of gravity.

Layered shallow water equations

1D and 2D solvers for the layered shallow water equations are also included.

Potentially useful contributions (what’s missing)

  • 2D mapped grid solvers (for a planar grid)

  • Transverse versions of rpn2_shallow_bathymetry_fwave.f90, rpn2_sw_aug.f90.

Demonstrations

[1]:
%matplotlib inline
import matplotlib.pyplot as plt
from clawpack import riemann
from clawpack.riemann import riemann_tools
import numpy as np
[2]:
h_l = 1.; h_r = 0.5;
u_l = 0.; u_r = 0.;
q_l = np.array([h_l,u_l]); q_r = np.array([h_r,u_r])
problem_data={'grav':1.0,'efix':False}

Roe

[3]:
states, speeds, reval = riemann_tools.riemann_solution(riemann.shallow_1D_py.shallow_roe_1D,q_l,q_r,
                                                       problem_data=problem_data)
riemann_tools.plot_phase(states)
riemann/../../../../../../../../Users/rjl/clawpack_src/clawpack_master/doc/doc/_build/html/v5.5.0/.doctrees/nbsphinx/riemann_Shallow_water_Riemann_solvers_9_0.png

HLL

[4]:
states, speeds, reval = riemann_tools.riemann_solution(riemann.shallow_1D_py.shallow_hll_1D,q_l,q_r,
                                                       problem_data=problem_data)
riemann_tools.plot_phase(states)
riemann/../../../../../../../../Users/rjl/clawpack_src/clawpack_master/doc/doc/_build/html/v5.5.0/.doctrees/nbsphinx/riemann_Shallow_water_Riemann_solvers_11_0.png