DispersiveShallowWater.jl API

DispersiveShallowWater.DispersiveShallowWaterModule
DispersiveShallowWater

DispersiveShallowWater.jl is a Julia package that implements structure-preserving numerical methods for dispersive shallow water models. It provides provably conservative, entropy-conserving, and well-balanced numerical schemes for some dispersive shallow water models.

The semidiscretizations are based on summation-by-parts (SBP) operators, which are implemented in SummationByPartsOperators.jl. To obtain fully discrete schemes, the time integration methods from OrdinaryDiffEq.jl are used to solve the resulting ordinary differential equations. Fully discrete entropy-conservative methods can be obtained by using the relaxation method provided by DispersiveShallowWater.jl.

See also: DispersiveShallowWater.jl

source

Equations

DispersiveShallowWater.BBMEquation1DType
BBMEquation1D(; gravity_constant, D = 1.0, eta0 = 0.0, split_form = true)

BBM (Benjamin–Bona–Mahony) equation in one spatial dimension. The equation is given by

\[\begin{aligned} \eta_t + \sqrt{gD}\eta_x + \frac{3}{2}\sqrt{\frac{g}{D}}\eta\eta_x - \frac{1}{6}D^2\eta_{xxt} &= 0. \end{aligned}\]

The unknown quantity of the BBM equation is the total water height $\eta$. The gravitational constant is denoted by g and the constant bottom topography (bathymetry) $b = \eta_0 - D$, where $\eta_0$ is the constant still-water surface and $D$ the still-water depth. The water height above the bathymetry is therefore given by $h = \eta - \eta_0 + D$. The BBM equation is only implemented for $\eta_0 = 0$.

The equations only support a flat bathymetry.

The BBM equation is first described in Benjamin, Bona, and Mahony (1972). The semidiscretization implemented here is developed in Ranocha, Mitsotakis, and Ketcheson (2020) for split_form = true and in Linders, Ranocha, and Birken (2023) for split_form = false. If split_form is true, a split form in the semidiscretization is used, which conserves

  • the total water mass (integral of $h$) as a linear invariant
  • a quadratic invariant (integral of $1/2\eta(\eta - 1/6D^2\eta_{xx})$ or for periodic boundary conditions equivalently $1/2(\eta^2 + 1/6D^2\eta_x^2)$), which is called here energy_total_modified (and entropy_modified) because it contains derivatives of the solution

for periodic boundary conditions. If split_form is false the semidiscretization conserves

  • the total water mass (integral of $h$) as a linear invariant
  • the Hamiltonian (integral of $1/4\sqrt{g/D}\eta^3 + 1/2\sqrt{gD}\eta^2$) (see hamiltonian)

for periodic boundary conditions.

  • Thomas B. Benjamin, Jerry L. Bona and John J. Mahony (1972) Model equations for long waves in nonlinear dispersive systems DOI: 10.1098/rsta.1972.0032
  • Hendrik Ranocha, Dimitrios Mitsotakis and David I. Ketcheson (2020) A Broad Class of Conservative Numerical Methods for Dispersive Wave Equations DOI: 10.4208/cicp.OA-2020-0119
  • Viktor Linders, Hendrik Ranocha and Philipp Birken (2023) Resolving entropy growth from iterative methods DOI: 10.1007/s10543-023-00992-w
source
DispersiveShallowWater.energy_total_modified!Method
energy_total_modified!(e, q_global, equations::BBMEquation1D, cache)

Return the modified total energy e of the primitive variables q_global for the BBMEquation1D. The energy_total_modified is a conserved quantity (for periodic boundary conditions).

It is given by

\[\frac{1}{2} \eta(\eta - \frac{1}{6}D^2\eta_{xx}).\]

q_global is a vector of the primitive variables at ALL nodes. cache needs to hold the SBP operators used by the solver.

See also energy_total_modified.

source
DispersiveShallowWater.hamiltonian!Method
hamiltonian!(H, q_global, equations::BBMEquation1D, cache)

Return the Hamiltonian H of the primitive variables q_global for the BBMEquation1D. The Hamiltonian is given by

\[\frac{1}{4}\sqrt{\frac{g}{D}}\eta^3 + \frac{1}{2}\sqrt{gD}\eta^2.\]

q_global is a vector of the primitive variables at ALL nodes.

See also hamiltonian.

source
DispersiveShallowWater.initial_condition_convergence_testMethod
initial_condition_convergence_test(x, t, equations::BBMEquation1D, mesh)

A travelling-wave solution used for convergence tests in a periodic domain, here generalized for dimensional variables.

See section 4.1.3 in (there is an error in paper: it should be sech^2 instead of cosh):

  • Hendrik Ranocha, Dimitrios Mitsotakis and David I. Ketcheson (2020) A Broad Class of Conservative Numerical Methods for Dispersive Wave Equations DOI: 10.4208/cicp.OA-2020-0119
source
DispersiveShallowWater.BBMBBMEquations1DType
BBMBBMEquations1D(; bathymetry_type = bathymetry_variable,
                  gravity_constant, eta0 = 0.0)

BBM-BBM (Benjamin–Bona–Mahony) system in one spatial dimension. The equations for flat bathymetry are given by

\[\begin{aligned} \eta_t + ((\eta + D)v)_x - \frac{1}{6}D^2\eta_{xxt} &= 0,\\ v_t + g\eta_x + \left(\frac{1}{2}v^2\right)_x - \frac{1}{6}D^2v_{xxt} &= 0. \end{aligned}\]

The unknown quantities of the BBM-BBM equations are the total water height $\eta$ and the velocity $v$. The gravitational constant is denoted by g and the constant bottom topography (bathymetry) $b = \eta_0 - D$. The water height above the bathymetry is therefore given by $h = \eta - \eta_0 + D$. The BBM-BBM equations are only implemented for $\eta_0 = 0$.

Two types of bathymetry_type are supported:

For the general case of variable vathymetry the BBM-BBM equations are

\[\begin{aligned} \eta_t + ((\eta + D)v)_x - \frac{1}{6}(D^2\eta_{xt})_x &= 0,\\ v_t + g\eta_x + \left(\frac{1}{2}v^2\right)_x - \frac{1}{6}(D^2v_t)_{xx} &= 0. \end{aligned}\]

One reference for the BBM-BBM system can be found in Bona et al. (1998). The semidiscretization implemented here was developed for flat bathymetry in Ranocha et al. (2020) and generalized for a variable bathymetry in Lampert and Ranocha (2024). It conserves

  • the total water mass (integral of $h$) as a linear invariant
  • the total velocity (integral of $v$) as a linear invariant for flat bathymetry
  • the total energy

for periodic boundary conditions (see Lampert, Ranocha). For reflecting boundary conditions, the semidiscretization conserves

  • the total water (integral of $h$) as a linear invariant
  • the total energy.

Additionally, it is well-balanced for the lake-at-rest stationary solution, see Lampert and Ranocha (2024).

  • Jerry L. Bona, Min Chen (1998) A Boussinesq system for two-way propagation of nonlinear dispersive waves DOI: 10.1016/S0167-2789(97)00249-2
  • Hendrik Ranocha, Dimitrios Mitsotakis, David I. Ketcheson (2020) A Broad Class of Conservative Numerical Methods for Dispersive Wave Equations DOI: 10.4208/cicp.OA-2020-0119
  • Joshua Lampert, Hendrik Ranocha (2024) Structure-Preserving Numerical Methods for Two Nonlinear Systems of Dispersive Wave Equations DOI: 10.48550/arXiv.2402.16669
source
DispersiveShallowWater.initial_condition_dingemansMethod
initial_condition_dingemans(x, t, equations::BBMBBMEquations1D, mesh)

The initial condition that uses the dispersion relation of the Euler equations to approximate waves generated by a wave maker as it is done by experiments of Dingemans. The topography is a trapezoidal.

Translation of water height

The initial condition for the water height is translated to be around 0, which is needed for the simulation because the BBMBBMEquations1D are only implemented for $\eta_0 = 0$.

References:

  • Magnus Svärd, Henrik Kalisch (2023) A novel energy-bounded Boussinesq model and a well-balanced and stable numerical discretization arXiv: 2302.09924
  • Maarten W. Dingemans (1994) Comparison of computations with Boussinesq-like models and laboratory measurements link
source
DispersiveShallowWater.AbstractShallowWaterEquationsType
AbstractShallowWaterEquations{NDIMS, NVARS}

An abstract supertype of all equation system that contain the classical shallow water equations as a subsystem, e.g., the BBMBBMEquations1D, the SvaerdKalischEquations1D, and the SerreGreenNaghdiEquations1D. In 1D, the shallow water equations with flat bathymetry are given by

\[\begin{aligned} h_t + (h v)_x &= 0,\\ h v_t + \frac{1}{2} g (h^2)_x + \frac{1}{2} h (v^2)_x &= 0, \end{aligned}\]

where $h$ is the waterheight, $v$ the velocity, and $g$ the gravity_constant.

source
DispersiveShallowWater.bathymetryMethod
bathymetry(q, equations)

Return the bathymetry of the primitive variables q for a given set of equations.

q is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations).

source
DispersiveShallowWater.cons2primMethod
cons2prim(u, equations)

Convert the conserved variables u to the primitive variables for a given set of equations. u is a vector type of the correct length nvariables(equations). Notice the function doesn't include any error checks for the purpose of efficiency, so please make sure your input is correct. The inverse conversion is performed by prim2cons.

source
DispersiveShallowWater.eachvariableMethod
eachvariable(equations::AbstractEquations)

Return an iterator over the indices that specify the location in relevant data structures for the variables in equations. In particular, not the variables themselves are returned.

source
DispersiveShallowWater.energy_totalMethod
energy_total(q, equations)

Return the total energy of the primitive variables q for a given set of equations. For all AbstractShallowWaterEquations, the total energy is given by the sum of the kinetic and potential energy of the shallow water subsystem, i.e.,

\[\frac{1}{2} h v^2 + \frac{1}{2} g \eta^2\]

in 1D, where $h$ is the waterheight, $\eta = h + b$ the waterheight_total, $v$ the velocity, and $g$ the gravity_constant.

q is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations).

source
DispersiveShallowWater.energy_total_modifiedMethod
energy_total_modified(q_global, equations, cache)

Return the modified total energy of the primitive variables q_global for the equations. This modified total energy is a conserved quantity and can contain additional terms compared to the usual energy_total. For example, for the SvaerdKalischEquations1D and the SerreGreenNaghdiEquations1D, it contains additional terms depending on the derivative of the velocity $v_x$ modeling non-hydrostatic contributions. For equations which are not AbstractShallowWaterEquations, the modified total energy does not have to be an extension of the usual energy_total and does not have to be related to a physical energy. However, it is still a conserved quantity for appropriate boundary conditions and may contain derivatives of the solution.

q_global is a vector of the primitive variables at ALL nodes. cache needs to hold the SBP operators used by the solver if non-hydrostatic terms are present.

Internally, this function allocates a vector for the output and calls DispersiveShallowWater.energy_total_modified!.

source
DispersiveShallowWater.get_nameMethod
get_name(equations::AbstractEquations)

Return the canonical, human-readable name for the given system of equations.

Examples

julia> DispersiveShallowWater.get_name(BBMBBMEquations1D(gravity_constant=1.0))
"BBMBBMEquations1D-BathymetryVariable"
source
DispersiveShallowWater.hamiltonianMethod
hamiltonian(q_global, equations, cache)

Return the Hamiltonian of the primitive variables q_global for the equations. The Hamiltonian is a conserved quantity and may contain derivatives of the solution.

q_global is a vector of the primitive variables at ALL nodes. cache needs to hold the SBP operators used by the solver if non-hydrostatic terms are present.

Internally, this function allocates a vector for the output and calls DispersiveShallowWater.hamiltonian!.

Note

This function is not necessarily implemented for all equations. See DispersiveShallowWater.hamiltonian! for further details of equations supporting it.

source
DispersiveShallowWater.hyperbolic_approximation_limitMethod
DispersiveShallowWater.hyperbolic_approximation_limit(equations)

If the equations are a hyperbolic approximation of another set of equations, return the equations of the limit system. Otherwise, return the input equations.

See also is_hyperbolic_appproximation and prim2phys.

Implementation details

This function is mostly used for some internal dispatch. For example, it allows to return a reduced set of variables from initial conditions for hyperbolic approximations.

source
DispersiveShallowWater.initial_condition_dingemansMethod
initial_condition_dingemans(x, t, equations::AbstractShallowWaterEquations, mesh)

The initial condition that uses the dispersion relation of the Euler equations to approximate waves generated by a wave maker as it is done by experiments of Dingemans. The topography is a trapezoidal. It is assumed that equations.eta0 = 0.8.

References:

  • Magnus Svärd, Henrik Kalisch (2023) A novel energy-bounded Boussinesq model and a well-balanced and stable numerical discretization arXiv: 2302.09924
  • Maarten W. Dingemans (1994) Comparison of computations with Boussinesq-like models and laboratory measurements link
source
DispersiveShallowWater.initial_condition_discontinuous_well_balancednessMethod
initial_condition_discontinuous_well_balancedness(x, t, equations::AbstractShallowWaterEquations, mesh)

Setup a truly discontinuous bottom topography function for this academic lake-at-rest test case of well-balancedness, i.e. eta is constant and v is zero everywhere. The error for this lake-at-rest test case ∫|η-η₀| should be around machine roundoff.

source
DispersiveShallowWater.is_hyperbolic_appproximationMethod
DispersiveShallowWater.is_hyperbolic_appproximation(equations)

Returns Val{true}() if the equations are a hyperbolic approximation of another set of equations and Val{false}() otherwise (default). For example, the HyperbolicSerreGreenNaghdiEquations1D are a hyperbolic approximation of the SerreGreenNaghdiEquations1D.

See also hyperbolic_approximation_limit and prim2phys.

Implementation details

This function is mostly used for some internal dispatch. For example, it allows to return a reduced set of variables from initial conditions for hyperbolic approximations.

source
DispersiveShallowWater.momentumMethod
momentum(q, equations)

Return the momentum/discharge of the primitive variables q for a given set of equations, i.e., the waterheight times the velocity.

q is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations).

source
DispersiveShallowWater.prim2consMethod
prim2cons(q, equations)

Convert the primitive variables q to the conserved variables for a given set of equations. q is a vector type of the correct length nvariables(equations). Notice the function doesn't include any error checks for the purpose of efficiency, so please make sure your input is correct. The inverse conversion is performed by cons2prim.

source
DispersiveShallowWater.prim2physMethod
prim2phys(q, equations)

Convert the primitive variables q to the physically meaningful variables for a given set of equations. By default, this is the same as prim2prim for most equations. However, some equations like the HyperbolicSerreGreenNaghdiEquations1D return a reduced set of variables since they are a hyperbolic approximation of another set of equations (in this case the SerreGreenNaghdiEquations1D).

See also is_hyperbolic_appproximation and hyperbolic_approximation_limit.

q is a vector type of the correct length nvariables(equations). Notice the function doesn't include any error checks for the purpose of efficiency, so please make sure your input is correct.

source
DispersiveShallowWater.velocityMethod
velocity(q, equations)

Return the velocity of the primitive variables q for a given set of equations.

q is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations).

source
DispersiveShallowWater.waterheightMethod
waterheight(q, equations)

Return the waterheight of the primitive variables q for a given set of equations, i.e., the waterheight $h$ above the bathymetry $b$.

q is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations).

See also waterheight_total, bathymetry.

source
DispersiveShallowWater.waterheight_totalMethod
waterheight_total(q, equations)

Return the total waterheight of the primitive variables q for a given set of equations, i.e., the waterheight $h$ plus the bathymetry $b$.

q is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations).

source
DispersiveShallowWater.HyperbolicSerreGreenNaghdiEquations1DType
HyperbolicSerreGreenNaghdiEquations1D(; bathymetry_type = bathymetry_mild_slope,
                                      gravity_constant,
                                      eta0 = 0.0,
                                      lambda)

Hyperbolic approximation of the Serre-Green-Naghdi system in one spatial dimension. The equations for flat bathymetry are given by

\[\begin{aligned} h_t + (h v)_x &= 0,\\ h v_t + \frac{1}{2} g (h^2)_x + \frac{1}{2} h (v^2)_x + \biggl( \frac{\lambda}{3} H (1 - H / h) \biggr)_x &= 0,\\ h w_t + h v w_x &= \lambda (1 - H / h),\\ H_t + H_x u &= w. \end{aligned}\]

The unknown quantities of the hyperbolized Serre-Green-Naghdi equations are the total water height $\eta = h + b$ and the velocity $v$. The gravitational constant is denoted by g and the bottom topography (bathymetry) $b = \eta_0 - D$. The water height above the bathymetry is therefore given by $h = \eta - \eta_0 + D$. The total water height is therefore given by $\eta = h + b$.

There are two additional variables $w \approx -h v_x$ and $H \approx h$ compared to the SerreGreenNaghdiEquations1D. In the original papers of Gavrilyuk et al., the variable $H$ is called $\eta$. Here, we use $\eta$ for the total water height and $H$ for auxiliary variable introduced in the hyperbolic approximation.

Initial conditions

The HyperbolicSerreGreenNaghdiEquations1D allow two options for specifying initial conditions:

  1. Returning the full set of variables q = (η, v, D, w, H)
  2. Returning a reduced set of variables q = (η, v, D) as required for the limit system SerreGreenNaghdiEquations1D of the hyperbolic approximation. The remaining variables w and H are initialized using the default initialization $w \approx -h v_x$ and $H \approx h$ using the derivative operator of the solver.

The relaxation parameter lambda ($\lambda$) introduced to obtain this hyperbolic approximation of the SerreGreenNaghdiEquations1D influences the stiffness of the system. For $\lambda \to \infty$, the hyperbolic Serre-Green-Naghdi equations converge (at least formally) to the original SerreGreenNaghdiEquations1D. However, the wave speeds of the hyperbolic system increase with increasing $\lambda$, so that explicit time integration methods become more expensive.

Two types of bathymetry_type are supported:

For the mild-slope approximation, the Serre-Green-Naghdi equations are

\[\begin{aligned} h_t + (h v)_x &= 0,\\ h v_t + \frac{1}{2} g (h^2)_x + \frac{1}{2} h (v^2)_x + \biggl( \frac{\lambda}{3} H (1 - H / h) \biggr)_x + \biggl( g h + \frac{\lambda}{2} (1 - H / h) \biggr) b_x &= 0,\\ h w_t + h v w_x &= \lambda (1 - H / h),\\ H_t + H_x u + \frac{3}{2} b_x v &= w. \end{aligned}\]

References for the hyperbolized Serre-Green-Naghdi system can be found in

  • Favrie and Gavrilyuk. A rapid numerical method for solving Serre-Green-Naghdi equations describing long free surface gravity waves DOI: 10.1088/1361-6544/aa712d
  • Busto, Dumbser, Escalante, Favrie, and Gavrilyuk. On High Order ADER Discontinuous Galerkin Schemes for First Order Hyperbolic Reformulations of Nonlinear Dispersive Systems DOI: 10.1007/s10915-021-01429-8

The semidiscretization implemented here conserves

  • the total water mass (integral of $h$) as a linear invariant
  • the total modified energy

for periodic boundary conditions (see Ranocha and Ricchiuto (2024)). Additionally, it is well-balanced for the lake-at-rest stationary solution, see

  • Hendrik Ranocha and Mario Ricchiuto (2024) Structure-preserving approximations of the Serre-Green-Naghdi equations in standard and hyperbolic form arXiv: 2408.02665
source
DispersiveShallowWater.energy_total_modified!Method
DispersiveShallowWater.energy_total_modified!(e, q_global, equations::HyperbolicSerreGreenNaghdiEquations1D, cache)

Return the modified total energy e of the primitive variables q_global for the HyperbolicSerreGreenNaghdiEquations1D. It contains additional terms compared to the usual energy_total modeling non-hydrostatic contributions. The energy_total_modified is a conserved quantity (for periodic boundary conditions).

For a bathymetry_mild_slope (and a bathymetry_flat), the total modified energy is given by

\[\frac{1}{2} g \eta^2 + \frac{1}{2} h v^2 + \frac{1}{6} h w^2 + \frac{\lambda}{6} h (1 - \eta / h)^2.\]

q_global is a vector of the primitive variables at ALL nodes.

See also energy_total_modified.

source
DispersiveShallowWater.initial_condition_solitonMethod
initial_condition_soliton(x, t, equations::HyperbolicSerreGreenNaghdiEquations1D, mesh)

A soliton solution of the SerreGreenNaghdiEquations1D used for convergence tests in a periodic domain. This is physically the same as initial_condition_convergence_test for the SerreGreenNaghdiEquations1D. Please note that this is not an exact solution of the HyperbolicSerreGreenNaghdiEquations1D (only in the limit of the relaxation parameter $\lambda \to \infty$).

See also initial_condition_convergence_test.

source
DispersiveShallowWater.SerreGreenNaghdiEquations1DType
SerreGreenNaghdiEquations1D(; bathymetry_type = bathymetry_variable,
                            gravity_constant, eta0 = 0.0)

Serre-Green-Naghdi system in one spatial dimension. The equations for flat bathymetry are given by

\[\begin{aligned} h_t + (h v)_x &= 0,\\ h v_t - \frac{1}{3} (h^3 v_{tx})_x + \frac{1}{2} g (h^2)_x + \frac{1}{2} h (v^2)_x + p_x &= 0,\\ p &= \frac{1}{3} h^3 v_{x}^2 - \frac{1}{3} h^3 v v_{xx}. \end{aligned}\]

The unknown quantities of the Serre-Green-Naghdi equations are the total water height $\eta = h + b$ and the velocity $v$. The gravitational constant is denoted by g and the bottom topography (bathymetry) $b = \eta_0 - D$. The water height above the bathymetry is therefore given by $h = \eta - \eta_0 + D$. The total water height is therefore given by $\eta = h + b$.

Three types of bathymetry_type are supported:

For the mild-slope approximation, the Serre-Green-Naghdi equations are

\[\begin{aligned} h_t + (h v)_x &= 0,\\ h v_t - \frac{1}{3} (h^3 v_{tx})_x + \frac{1}{2} (h^2 b_x u_t)_x - \frac{1}{2} h^2 b_x u_{tx} + \frac{3}{4} h b_x^2 u_t + \frac{1}{2} g (h^2)_x + g h b_x + \frac{1}{2} h (v^2)_x + p_x + \frac{3}{2} \frac{p}{h} b_x &= 0,\\ p &= \frac{1}{3} h^3 v_{x}^2 - \frac{1}{3} h^3 v v_{xx} + \frac{1}{2} h^2 v (b_x v)_x. \end{aligned}\]

For the general case of variable vathymetry without mild-slope approximation, the Serre-Green-Naghdi equations are

\[\begin{aligned} h_t + (h v)_x &= 0,\\ h v_t - \frac{1}{3} (h^3 v_{tx})_x + \frac{1}{2} (h^2 b_x u_t)_x - \frac{1}{2} h^2 b_x u_{tx} + h b_x^2 u_t + \frac{1}{2} g (h^2)_x + g h b_x + \frac{1}{2} h (v^2)_x + p_x + \frac{3}{2} \frac{p}{h} b_x + \psi b_x &= 0,\\ p &= \frac{1}{3} h^3 v_{x}^2 - \frac{1}{3} h^3 v v_{xx} + \frac{1}{2} h^2 v (b_x v)_x,\\ \psi &= \frac{1}{4} h v (b_x v)_x. \end{aligned}\]

References for the Serre-Green-Naghdi system can be found in

The semidiscretization implemented here conserves

  • the total water mass (integral of $h$) as a linear invariant
  • the total momentum (integral of $h v$) as a nonlinear invariant if the bathymetry is constant
  • the total modified energy

for periodic boundary conditions (see Ranocha and Ricchiuto (2024)). Additionally, it is well-balanced for the lake-at-rest stationary solution, see

  • Hendrik Ranocha and Mario Ricchiuto (2024) Structure-preserving approximations of the Serre-Green-Naghdi equations in standard and hyperbolic form arXiv: 2408.02665
source
DispersiveShallowWater.energy_total_modified!Method
DispersiveShallowWater.energy_total_modified!(e, q_global, equations::SerreGreenNaghdiEquations1D, cache)

Return the modified total energy e of the primitive variables q_global for the SerreGreenNaghdiEquations1D. It contains an additional term containing a derivative compared to the usual energy_total modeling non-hydrostatic contributions. The energy_total_modified is a conserved quantity (for periodic boundary conditions).

For a bathymetry_flat the total modified energy is given by

\[\frac{1}{2} g \eta^2 + \frac{1}{2} h v^2 + \frac{1}{6} h^3 v_x^2.\]

For a bathymetry_mild_slope the total modified energy is given by

\[\frac{1}{2} g \eta^2 + \frac{1}{2} h v^2 + \frac{1}{6} h (-h v_x + 1.5 v b_x)^2.\]

For a bathymetry_variable the total modified energy has the additional term

\[+ \frac{1}{8} h (v b_x)^2.\]

q_global is a vector of the primitive variables at ALL nodes. cache needs to hold the SBP operators used by the solver.

See also energy_total_modified.

source
DispersiveShallowWater.SvaerdKalischEquations1DType
SvaerdKalischEquations1D(; bathymetry_type = bathymetry_variable,
                         gravity_constant, eta0 = 0.0, alpha = 0.0,
                         beta = 0.2308939393939394, gamma = 0.04034343434343434)

Dispersive system by Svärd and Kalisch (2023) in one spatial dimension. The equations for variable bathymetry are given in conservative variables by

\[\begin{aligned} h_t + (hv)_x &= (\hat\alpha(\hat\alpha(h + b)_x)_x)_x,\\ (hv)_t + (hv^2)_x + gh(h + b)_x &= (\hat\alpha v(\hat\alpha(h + b)_x)_x)_x + (\hat\beta v_x)_{xt} + \frac{1}{2}(\hat\gamma v_x)_{xx} + \frac{1}{2}(\hat\gamma v_{xx})_x, \end{aligned}\]

where $\hat\alpha^2 = \alpha\sqrt{gD}D^2$, $\hat\beta = \beta D^3$, $\hat\gamma = \gamma\sqrt{gD}D^3$. The coefficients $\alpha$, $\beta$ and $\gamma$ are provided in dimensionless form and $D = \eta_0 - b$ is the still-water depth and eta0 is the still-water surface (lake-at-rest). The equations can be rewritten in primitive variables as

\[\begin{aligned} \eta_t + ((\eta + D)v)_x &= (\hat\alpha(\hat\alpha\eta_x)_x)_x,\\ v_t(\eta + D) - v((\eta + D)v)_x + ((\eta + D)v^2)_x + g(\eta + D)\eta_x &= (\hat\alpha v(\hat\alpha\eta_x)_x)_x - v(\hat\alpha(\hat\alpha\eta_x)_x)_x + (\hat\beta v_x)_{xt} + \frac{1}{2}(\hat\gamma v_x)_{xx} + \frac{1}{2}(\hat\gamma v_{xx})_x. \end{aligned}\]

The unknown quantities of the Svärd-Kalisch equations are the total water height $\eta$ and the velocity $v$. The gravitational constant is denoted by g and the bottom topography (bathymetry) $b = \eta_0 - D$. The water height above the bathymetry is therefore given by $h = \eta - \eta_0 + D$.

Currently, the equations only support a general variable bathymetry, see bathymetry_variable.

SvärdKalischEquations1D is an alias for SvaerdKalischEquations1D.

The equations by Svärd and Kalisch are presented and analyzed in Svärd and Kalisch (2023). The semidiscretization implemented here conserves

  • the total water mass (integral of $h$) as a linear invariant
  • the total momentum (integral of $h v$) as a nonlinear invariant for flat bathymetry
  • the total modified energy

for periodic boundary conditions (see Lampert, Ranocha). Additionally, it is well-balanced for the lake-at-rest stationary solution, see Lampert and Ranocha (2024).

  • Magnus Svärd, Henrik Kalisch (2023) A novel energy-bounded Boussinesq model and a well-balanced and stable numerical discretization arXiv: 2302.09924
  • Joshua Lampert, Hendrik Ranocha (2024) Structure-Preserving Numerical Methods for Two Nonlinear Systems of Dispersive Wave Equations DOI: 10.48550/arXiv.2402.16669
source
DispersiveShallowWater.energy_total_modified!Method
DispersiveShallowWater.energy_total_modified!(e, q_global, equations::SvaerdKalischEquations1D, cache)

Return the modified total energy e of the primitive variables q_global for the SvaerdKalischEquations1D. It contains an additional term containing a derivative compared to the usual energy_total modeling non-hydrostatic contributions. The energy_total_modified is a conserved quantity (for periodic boundary conditions).

It is given by

\[\frac{1}{2} g \eta^2 + \frac{1}{2} h v^2 + \frac{1}{2} \hat\beta v_x^2.\]

q_global is a vector of the primitive variables at ALL nodes. cache needs to hold the SBP operators used by the solver.

See also energy_total_modified.

source

Mesh

Boundary conditions

Solver

DispersiveShallowWater.SolverMethod
Solver(mesh, accuracy_order)

Create a solver, where the summation-by-parts (SBP) operators are of order accuracy_order and associated to the mesh.

source
DispersiveShallowWater.SolverMethod
Solver(D1, D2)

Create a solver, where D1 is an AbstractDerivativeOperator from SummationByPartsOperators.jl of first derivative_order and D2 is an AbstractDerivativeOperator of second derivative_order or an AbstractMatrix. It can also be nothing if no second derivative is used by the discretization. Both summation-by-parts operators should be associated with the same grid.

source

Semidiscretization

DispersiveShallowWater.SemidiscretizationMethod
Semidiscretization(mesh, equations, initial_condition, solver;
                   source_terms=nothing,
                   boundary_conditions=boundary_condition_periodic,
                   RealT=real(solver),
                   uEltype=RealT,
                   initial_cache=(tmp1 = Array{RealT}(undef, nnodes(mesh)),))

Construct a semidiscretization of a PDE.

source

Callbacks

DispersiveShallowWater.AnalysisCallbackType
AnalysisCallback(semi; interval=0,
                       extra_analysis_errors=Symbol[],
                       extra_analysis_integrals=(),
                       io=stdout)

Analyze a numerical solution every interval time steps. The L2- and the L∞-norm for each component are computed by default. Additional errors can be computed, e.g. by passing extra_analysis_errors = (:conservation_error,).

Further scalar functions func in extra_analysis_integrals are applied to the numerical solution and integrated over the computational domain. Some examples for this are entropy, and energy_total. You can also write your own function with the same signature as the examples listed above and pass it via extra_analysis_integrals. The computed errors and intergrals are saved for each timestep and can be obtained by calling errors and integrals.

During the Simulation, the AnalysisCallback will print information to io.

source
DispersiveShallowWater.errorsMethod
errors(analysis_callback)

Return the computed errors for each timestep as a named tuple. The shape of each entry is (nvariables, ntimesteps).

source
DispersiveShallowWater.RelaxationCallbackType
RelaxationCallback(invariant)

Use a relaxation method in time in order to exactly preserve the (nonlinear) invariant for a conservative semidiscretization. A possible choice for invariant is invariant = entropy.

Reference

  • Hendrik Ranocha, Mohammed Sayyari, Lisandro Dalcin, Matteo Parsani, David I. Ketcheson (2020) Relaxation Runge–Kutta Methods: Fully-Discrete Explicit Entropy-Stable Schemes for the Compressible Euler and Navier–Stokes Equations DOI: 10.1137/19M1263480
source
DispersiveShallowWater.SummaryCallbackType
SummaryCallback(io::IO = stdout)

Create and return a callback that resets the timer at the beginning of a simulation and prints the timer values at the end of the simulation.

source

Utilities

DispersiveShallowWater.convergence_testMethod
convergence_test([mod::Module=Main,] example::AbstractString, iterations; io::IO = stdout, kwargs...)
convergence_test([mod::Module=Main,] example::AbstractString, Ns::AbstractVector; io::IO = stdout, kwargs...)

Run multiple simulations using the setup given in example and compute the experimental order of convergence (EOC) in the $L^2$ and $L^\infty$ norm. If iterations is passed as integer, in each iteration, the resolution of the respective mesh will be doubled. If Ns is passed as vector, the simulations will be run for each value of Ns. Additional keyword arguments kwargs... and the optional module mod are passed directly to trixi_include.

Adjusted from Trixi.jl.

source
DispersiveShallowWater.examples_dirMethod
examples_dir()

Return the directory where the example files provided with DispersiveShallowWater.jl are located. If DispersiveShallowWater is installed as a regular package (with ]add DispersiveShallowWater), these files are read-only and should not be modified. To find out which files are available, use, e.g., readdir.

Copied from Trixi.jl.

Examples

readdir(examples_dir())
source
DispersiveShallowWater.@autoinfiltrateMacro
@autoinfiltrate
@autoinfiltrate condition::Bool

Invoke the @infiltrate macro of the package Infiltrator.jl to create a breakpoint for ad-hoc interactive debugging in the REPL. If the optional argument condition is given, the breakpoint is only enabled if condition evaluates to true.

As opposed to using Infiltrator.@infiltrate directly, this macro does not require Infiltrator.jl to be added as a dependency to DispersiveShallowWater.jl. As a bonus, the macro will also attempt to load the Infiltrator module if it has not yet been loaded manually.

Note: For this macro to work, the Infiltrator.jl package needs to be installed in your current Julia environment stack.

See also: Infiltrator.jl

Internal use only

Please note that this macro is intended for internal use only. It is not part of the public API of DispersiveShallowWater.jl, and it thus can altered (or be removed) at any time without it being considered a breaking change.

source