DispersiveShallowWater.jl API
DispersiveShallowWater.DispersiveShallowWater
— ModuleDispersiveShallowWater
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
Model specific equations
DispersiveShallowWater.energy_total_modified!
— Methodenergy_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
.
DispersiveShallowWater.hamiltonian!
— Methodhamiltonian!(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
.
DispersiveShallowWater.initial_condition_convergence_test
— Methodinitial_condition_convergence_test(x, t, equations::BBMEquation1D, mesh)
A traveling-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
DispersiveShallowWater.initial_condition_manufactured
— Methodinitial_condition_manufactured(x, t, equations::BBMEquation1D, mesh)
A smooth manufactured solution in combination with source_terms_manufactured
.
DispersiveShallowWater.source_terms_manufactured
— Methodsource_terms_manufactured(q, x, t, equations::BBMEquation1D)
A smooth manufactured solution in combination with initial_condition_manufactured
.
DispersiveShallowWater.initial_condition_convergence_test
— Methodinitial_condition_convergence_test(x, t, equations::BBMBBMEquations1D, mesh)
A traveling-wave solution used for convergence tests in a periodic domain. The bathymetry is constant.
For details see Example 5 in Section 3 from (here adapted for dimensional equations):
- Min Chen (1997) Exact Traveling-Wave Solutions to Bidirectional Wave Equations DOI: 10.1023/A:1026667903256
DispersiveShallowWater.initial_condition_dingemans
— Methodinitial_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.
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
DispersiveShallowWater.initial_condition_manufactured
— Methodinitial_condition_manufactured(x, t, equations::BBMBBMEquations1D, mesh)
A smooth manufactured solution in combination with source_terms_manufactured
.
DispersiveShallowWater.initial_condition_manufactured_reflecting
— Methodinitial_condition_manufactured_reflecting(x, t, equations::BBMBBMEquations1D, mesh)
A smooth manufactured solution for reflecting boundary conditions in combination with source_terms_manufactured_reflecting
.
DispersiveShallowWater.source_terms_manufactured
— Methodsource_terms_manufactured(q, x, t, equations::BBMBBMEquations1D)
A smooth manufactured solution in combination with initial_condition_manufactured
.
DispersiveShallowWater.source_terms_manufactured_reflecting
— Methodsource_terms_manufactured_reflecting(q, x, t, equations::BBMBBMEquations1D)
A smooth manufactured solution for reflecting boundary conditions in combination with initial_condition_manufactured_reflecting
.
DispersiveShallowWater.bathymetry_flat
— Constantbathymetry_flat = DispersiveShallowWater.BathymetryFlat()
A singleton struct indicating a flat bathymetry.
See also bathymetry_mild_slope
and bathymetry_variable
.
DispersiveShallowWater.bathymetry_mild_slope
— Constantbathymetry_mild_slope = DispersiveShallowWater.BathymetryMildSlope()
A singleton struct indicating a variable bathymetry with mild-slope approximation. Typically, this means that some terms like $b_x^2$ are neglected.
See also bathymetry_flat
and bathymetry_variable
.
DispersiveShallowWater.bathymetry_variable
— Constantbathymetry_variable = DispersiveShallowWater.BathymetryVariable()
A singleton struct indicating a variable bathymetry (without mild-slope approximation).
See also bathymetry_flat
and bathymetry_mild_slope
.
DispersiveShallowWater.AbstractEquations
— TypeAbstractEquations{NDIMS, NVARS}
An abstract supertype of specific equations such as the BBM-BBM equations. The type parameters encode the number of spatial dimensions (NDIMS
) and the number of primary variables (NVARS
) of the physics model.
See also AbstractShallowWaterEquations
.
DispersiveShallowWater.AbstractShallowWaterEquations
— TypeAbstractShallowWaterEquations{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
.
DispersiveShallowWater.bathymetry
— Methodbathymetry(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)
.
DispersiveShallowWater.cons2prim
— Methodcons2prim(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
.
DispersiveShallowWater.default_analysis_errors
— Methoddefault_analysis_errors(equations)
Default analysis errors used by the AnalysisCallback
.
DispersiveShallowWater.default_analysis_integrals
— Methoddefault_analysis_integrals(equations)
Default analysis integrals used by the AnalysisCallback
.
DispersiveShallowWater.discharge
— Methoddischarge(q, equations)
See momentum
.
DispersiveShallowWater.eachvariable
— Methodeachvariable(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.
DispersiveShallowWater.energy_total
— Methodenergy_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
.
q
is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations)
.
DispersiveShallowWater.energy_total_modified!
— Methodenergy_total_modified!(e, q_global, equations, cache)
In-place version of energy_total_modified
.
DispersiveShallowWater.energy_total_modified
— Methodenergy_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!
.
DispersiveShallowWater.entropy
— Methodentropy(q, equations)
Return the entropy of the primitive variables q
for a given set of equations
. For all AbstractShallowWaterEquations
, the entropy
is just the energy_total
.
q
is a vector of the primitive variables at a single node, i.e., a vector of the correct length nvariables(equations)
.
DispersiveShallowWater.entropy_modified!
— Methodentropy_modified!(e, q_global, equations, cache)
In-place version of entropy_modified
.
DispersiveShallowWater.entropy_modified
— Methodentropy_modified(q_global, equations, cache)
Alias for energy_total_modified
.
DispersiveShallowWater.get_name
— Methodget_name(equations::AbstractEquations)
Return the canonical, human-readable name for the given system of equations.
Examples
julia> DispersiveShallowWater.get_name(BBMBBMEquations1D(gravity=1.0))
"BBMBBMEquations1D-BathymetryVariable"
DispersiveShallowWater.gravity
— Methodgravity(equations)
Return the gravitational acceleration $g$ for a given set of equations
.
DispersiveShallowWater.hamiltonian
— Methodhamiltonian(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!
.
This function is not necessarily implemented for all equations. See DispersiveShallowWater.hamiltonian!
for further details of equations supporting it.
DispersiveShallowWater.hyperbolic_approximation_limit
— MethodDispersiveShallowWater.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
.
DispersiveShallowWater.initial_condition_dingemans
— Methodinitial_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
DispersiveShallowWater.initial_condition_discontinuous_well_balancedness
— Methodinitial_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.
DispersiveShallowWater.is_hyperbolic_appproximation
— MethodDispersiveShallowWater.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
.
DispersiveShallowWater.lake_at_rest_error
— Methodlake_at_rest_error(q, equations)
Calculate the error for the "lake-at-rest" test case where the waterheight_total
$\eta = h + b$ should be a constant value over time (given by the value $\eta_0$ passed to the equations
when constructing them).
DispersiveShallowWater.momentum
— Methodmomentum(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)
.
DispersiveShallowWater.prim2cons
— Methodprim2cons(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
.
DispersiveShallowWater.prim2phys
— Methodprim2phys(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.
DispersiveShallowWater.prim2prim
— Methodprim2prim(q, equations)
Return the primitive variables q
. While this function is as trivial as identity
, it is also as useful.
DispersiveShallowWater.still_water_surface
— Methodstill_water_surface(q, equations)
Return the still water surface $\eta_0$ (lake at rest) for a given set of equations
.
DispersiveShallowWater.varnames
— Functionvarnames(conversion_function, equations)
Return the list of variable names when applying conversion_function
to the primitive variables associated to equations
. Common choices of the conversion_function
are prim2prim
, prim2cons
, and prim2phys
.
DispersiveShallowWater.velocity
— Methodvelocity(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)
.
DispersiveShallowWater.waterheight
— Methodwaterheight(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
.
DispersiveShallowWater.waterheight_total
— Methodwaterheight_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)
.
DispersiveShallowWater.energy_total_modified!
— MethodDispersiveShallowWater.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
.
DispersiveShallowWater.initial_condition_manufactured
— Methodinitial_condition_manufactured(x, t, equations::HyperbolicSerreGreenNaghdiEquations1D, mesh)
A smooth manufactured solution in combination with source_terms_manufactured
, see
- Hendrik Ranocha and Mario Ricchiuto (2024) Structure-preserving approximations of the Serre-Green-Naghdi equations in standard and hyperbolic form arXiv: 2408.02665
DispersiveShallowWater.initial_condition_manufactured_reflecting
— Methodinitial_condition_manufactured_reflecting(x, t, equations::HyperbolicSerreGreenNaghdiEquations1D, mesh)
A smooth manufactured solution for reflecting boundary conditions in combination with source_terms_manufactured_reflecting
.
DispersiveShallowWater.initial_condition_soliton
— Methodinitial_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
.
DispersiveShallowWater.prim2phys
— Methodprim2phys(q, equations::HyperbolicSerreGreenNaghdiEquations1D)
Return the physical variables $\eta, v, D$ used also by the SerreGreenNaghdiEquations1D
from the main variables q
for the HyperbolicSerreGreenNaghdiEquations1D
.
DispersiveShallowWater.source_terms_manufactured
— Methodsource_terms_manufactured(q, x, t, equations::HyperbolicSerreGreenNaghdiEquations1D)
A smooth manufactured solution in combination with initial_condition_manufactured
.
DispersiveShallowWater.source_terms_manufactured_reflecting
— Methodsource_terms_manufactured_reflecting(q, x, t, equations::HyperbolicSerreGreenNaghdiEquations1D)
A smooth manufactured solution in combination with initial_condition_manufactured_reflecting
.
calculated as shown in: https://github.com/NumericalMathematics/DispersiveShallowWater.jl/pull/228#issuecomment-3123350726
DispersiveShallowWater.energy_total
— Methodenergy_total(q, equations::KdVEquation1D)
Return the total energy of the primitive variables q
for the KdVEquation1D
. For the KdV equation, the total energy consists only of the potential energy, given by
\[\frac{1}{2} g \eta^2\]
where $\eta$ is the waterheight_total
and $g$ is the gravity
.
q
is a vector of the primitive variables at a single node, i.e., a vector of length 1 in this case.
DispersiveShallowWater.entropy
— Methodentropy(q, equations::KdVEquation1D)
Return the entropy of the primitive variables q
for the KdVEquation1D
. For the KdV equation, the entropy
is the same as the energy_total
.
q
is a vector of the primitive variables at a single node, i.e., a vector of length 1 in this case.
DispersiveShallowWater.initial_condition_convergence_test
— Methodinitial_condition_convergence_test(x, t, equations::KdVEquation1D, mesh)
A soliton solution used for convergence tests in a periodic domain. Same as initial_condition_soliton
for the KdVEquation1D
.
DispersiveShallowWater.initial_condition_manufactured
— Methodinitial_condition_manufactured(x, t, equations::KdVEquation1D, mesh)
A smooth manufactured solution in combination with source_terms_manufactured
.
DispersiveShallowWater.initial_condition_soliton
— Methodinitial_condition_soliton(x, t, equations::KdVEquation1D, mesh)
A classical soliton solution of the KdV equation in dimensional variables. This can be used for convergence tests in a periodic domain, see initial_condition_convergence_test
.
DispersiveShallowWater.nondim2prim
— Methodnondim2prim(u, equations::KdVEquation1D)
Convert the non-dimensional variable u
to the primitive/physical variable eta
(total water height) for the KdVEquation1D
.
The transformation is given by:
\[\eta = D(u - \frac{2}{3})\]
where D
is the still-water depth.
This conversion is only valid for equations with specific parameter values:
gravity = 4/27
D = 3.0
These values ensure the dimensional KdV equation matches the standard non-dimensional form u_t + u u_x + u_{xxx} = 0
.
This function allows converting solutions from the standard non-dimensional KdV form commonly found in literature to the dimensional form implemented in DispersiveShallowWater.jl.
See also prim2nondim
.
DispersiveShallowWater.prim2nondim
— Methodprim2nondim(eta, equations::KdVEquation1D)
Convert the primitive/physical variable eta
(total water height) to the non-dimensional variable u
for the KdVEquation1D
.
The transformation is given by:
\[u = \frac{\eta}{D} + \frac{2}{3}\]
where D
is the still-water depth.
This conversion is only valid for equations with specific parameter values:
gravity = 4/27
D = 3.0
These values ensure the dimensional KdV equation matches the standard non-dimensional form u_t + u u_x + u_{xxx} = 0
.
This function allows converting solutions from the dimensional form implemented in DispersiveShallowWater.jl to the standard non-dimensional KdV form commonly found in literature, enabling comparison with theoretical results and other implementations.
See also nondim2prim
.
DispersiveShallowWater.source_terms_manufactured
— Methodsource_terms_manufactured(q, x, t, equations::KdVEquation1D)
A smooth manufactured solution in combination with initial_condition_manufactured
.
How it was calculated, is described in: https://github.com/NumericalMathematics/DispersiveShallowWater.jl/pull/198#discussion_r2090805751
DispersiveShallowWater.varnames
— Methodvarnames(::typeof(prim2nondim), equations::KdVEquation1D)
Return variable names ("u",)
for non-dimensional KdV variables when plotting with conversion = prim2nondim
.
See prim2nondim
, varnames
.
DispersiveShallowWater.energy_total_modified!
— MethodDispersiveShallowWater.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
.
DispersiveShallowWater.initial_condition_convergence_test
— Methodinitial_condition_convergence_test(x, t, equations::SerreGreenNaghdiEquations1D, mesh)
A soliton solution used for convergence tests in a periodic domain. Same as initial_condition_soliton
for the SerreGreenNaghdiEquations1D
.
DispersiveShallowWater.initial_condition_manufactured
— Methodinitial_condition_manufactured(x, t, equations::SerreGreenNaghdiEquations1D, mesh)
A smooth manufactured solution in combination with source_terms_manufactured
.
DispersiveShallowWater.initial_condition_manufactured_reflecting
— Methodinitial_condition_manufactured_reflecting(x, t, equations::SerreGreenNaghdiEquations1D, mesh)
A smooth manufactured solution for reflecting boundary conditions in combination with source_terms_manufactured_reflecting
.
References for flat bathymetry
- D. Mitsotakis, C. Synolakis, and M. McGuinness (2016) A modified Galerkin/finite element method for the numerical solution of the Serre-Green-Naghdi system. DOI: 10.1002/fld.4293
DispersiveShallowWater.initial_condition_soliton
— Methodinitial_condition_soliton(x, t, equations::SerreGreenNaghdiEquations1D, mesh)
A classical soliton solution of the Serre-Green-Naghdi equations. This can be used for convergence tests in a periodic domain, see initial_condition_convergence_test
.
DispersiveShallowWater.source_terms_manufactured
— Methodsource_terms_manufactured(q, x, t, equations::SerreGreenNaghdiEquations1D)
A smooth manufactured solution in combination with initial_condition_manufactured
.
DispersiveShallowWater.source_terms_manufactured_reflecting
— Methodsource_terms_manufactured_reflecting(q, x, t, equations::SerreGreenNaghdiEquations1D)
A smooth manufactured solution for reflecting boundary conditions in combination with initial_condition_manufactured_reflecting
.
DispersiveShallowWater.energy_total_modified!
— MethodDispersiveShallowWater.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
.
DispersiveShallowWater.initial_condition_manufactured
— Methodinitial_condition_manufactured(x, t, equations::SvaerdKalischEquations1D, mesh)
A smooth manufactured solution in combination with source_terms_manufactured
.
DispersiveShallowWater.initial_condition_manufactured_reflecting
— Methodinitial_condition_manufactured_reflecting(x, t, equations::SvaerdKalischEquations1D, mesh)
A smooth manufactured solution for reflecting boundary conditions in combination with source_terms_manufactured_reflecting
.
DispersiveShallowWater.source_terms_manufactured
— Methodsource_terms_manufactured(q, x, t, equations::SvaerdKalischEquations1D)
A smooth manufactured solution in combination with initial_condition_manufactured
.
DispersiveShallowWater.source_terms_manufactured_reflecting
— Methodsource_terms_manufactured_reflecting(q, x, t, equations::SvaerdKalischEquations1D)
A smooth manufactured solution for reflecting boundary conditions in combination with initial_condition_manufactured_reflecting
.
Linear dispersion relations
DispersiveShallowWater.EulerEquations1D
— TypeEulerEquations1D(; gravity, eta0 = 0.0)
A struct representing the 1D Euler equations with a given gravitational acceleration gravity
and a still-water surface eta0
.
In DispersiveShallowWater.jl, the Euler equations are only used for computing the full linear dispersion relation
\[\omega(k) = \sqrt{g k \tanh(h_0 k)},\]
see LinearDispersionRelation
and wave_speed
. The EulerEquations1D
cannot be solved as a system of equations.
DispersiveShallowWater.LinearDispersionRelation
— TypeLinearDispersionRelation(ref_height)
A struct representing a linear dispersion relation $\omega(k)$ of an equation. The reference water height $h_0$ is given by ref_height
. A dispersion relation can be called as disp_rel(equations, k)
to compute the wave frequency $\omega(k)$ for a given wavenumber k
and a set of equations.
See also wave_speed
for computing the wave speed $c = \omega(k) / k$ given a linear dispersion relation.
DispersiveShallowWater.wave_speed
— Methodwave_speed(disp_rel, equations, k; normalize = false)
Compute the wave speed $c$ for a given wavenumber $k$ using the LinearDispersionRelation
disp_rel
of the equations
. The wave speed is given by $c = \omega(k) / k$. If normalize
is true
, the wave speed is normalized by the shallow water wave speed $\sqrt{g h_0}$, where $g$ is the gravitational acceleration of the equations
and $h_0$ is the ref_height
of the dispersion relation disp_rel
.
See also LinearDispersionRelation
.
Mesh
DispersiveShallowWater.Mesh1D
— TypeMesh1D
Struct that holds the information for a simple homogeneous one-dimensional mesh.
DispersiveShallowWater.Mesh1D
— MethodMesh1D(xmin, xmax, N)
Create a simple homogeneous one-dimensional mesh from xmin
to xmax
with N
nodes.
Boundary conditions
DispersiveShallowWater.boundary_condition_periodic
— Constantboundary_condition_periodic = DispersiveShallowWater.BoundaryConditionPeriodic()
A singleton struct indicating periodic boundary conditions.
DispersiveShallowWater.boundary_condition_reflecting
— Constantboundary_condition_reflecting = DispersiveShallowWater.BoundaryConditionReflecting()
A singleton struct indicating reflecting boundary conditions.
Solver
DispersiveShallowWater.AbstractSolver
— TypeAbstractSolver
An abstract supertype of specific solvers.
DispersiveShallowWater.Solver
— TypeSolver
A struct
that holds the summation-by-parts (SBP) operators that are used for the spatial discretization.
DispersiveShallowWater.Solver
— MethodSolver(mesh, accuracy_order)
Create a solver, where the three summation-by-parts (SBP) first-, second-, and third-derivative operators are of accuracy order accuracy_order
and associated to the mesh
.
DispersiveShallowWater.Solver
— MethodSolver(D1, D2 = nothing, D3 = nothing)
Create a solver, where D1
is an AbstractDerivativeOperator
from SummationByPartsOperators.jl of first derivative_order
, D2
is an AbstractDerivativeOperator
of second derivative_order
or an AbstractMatrix
, and D3
is an AbstractDerivativeOperator
of third derivative_order
or an AbstractMatrix
. D2
and D3
can also be nothing
if that derivative is not used by the discretization. All given summation-by-parts operators should be associated with the same grid.
Semidiscretization
DispersiveShallowWater.Semidiscretization
— TypeSemidiscretization
A struct
containing everything needed to describe a spatial semidiscretization of an equation.
DispersiveShallowWater.Semidiscretization
— MethodSemidiscretization(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)),
tmp_partitioned = allocate_coefficients(mesh, equations, solver)))
Construct a semidiscretization of a PDE.
DispersiveShallowWater.check_solver
— Methodcheck_solver(equations, solver, boundary_conditions)
Check that the solver
is compatible with the given equations
and boundary_conditions
. The default implementation performs no checks. Specific equation types can override this method to validate that required derivative operators are present (e.g., some equations require D2
or D3
to be non-nothing
).
Throws an ArgumentError
if the solver is incompatible.
DispersiveShallowWater.jacobian
— MethodDispersiveShallowWater.jacobian(semi::Semidiscretization;
t = 0.0,
q0 = compute_coefficients(semi.initial_condition, t, semi))
Use the right-hand side operator of the semidiscretization semi
and forward mode automatic differentiation to compute the Jacobian J
of the semidiscretization semi
at the state q0
.
PolynomialBases.grid
— Methodgrid(semi)
Get the grid of a semidiscretization.
SummationByPartsOperators.semidiscretize
— Methodsemidiscretize(semi::Semidiscretization, tspan)
Wrap the semidiscretization semi
as an ODE problem in the time interval tspan
that can be passed to solve
from the SciML ecosystem.
Callbacks
DispersiveShallowWater.AnalysisCallback
— TypeAnalysisCallback(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
.
DispersiveShallowWater.errors
— Methoderrors(analysis_callback)
Return the computed errors for each timestep as a named tuple. The shape of each entry is (nvariables, ntimesteps).
DispersiveShallowWater.integrals
— Methodintegrals(analysis_callback)
Return the computed integrals for each timestep as a named tuple.
DispersiveShallowWater.tstops
— Methodtstops(analysis_callback)
Return the time values that correspond to the saved values of the errors
and integrals
.
DispersiveShallowWater.RelaxationCallback
— TypeRelaxationCallback(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
DispersiveShallowWater.SummaryCallback
— TypeSummaryCallback(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.
Utilities
Experimental data
DispersiveShallowWater.data_dingemans
— Methoddata_dingemans()
Load the experimental data for the Dingemans experiment. Returns the time values, x-coordinates of the six wave gauges, and experimental data, which is a matrix of wave heights at each gauge location (columns) over time (rows).
- Dingemans (1994) Comparison of computations with Boussinesq-like models and laboratory measurements URL: https://resolver.tudelft.nl/uuid:c2091d53-f455-48af-a84b-ac86680455e9
- Dingemans (1997): Water Wave Propagation Over Uneven Bottoms (In 2 Parts). DOI: 10.1142/1241
Utilities from DispersiveShallowWater.jl
DispersiveShallowWater.convergence_test
— Methodconvergence_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.
DispersiveShallowWater.data_dir
— Methoddata_dir()
Return the directory where the data files provided with DispersiveShallowWater.jl are located. If DispersiveShallowWater.jl 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
.
readdir(data_dir())
DispersiveShallowWater.default_example
— Methoddefault_example()
Return the path to an example that can be used to quickly see DispersiveShallowWater.jl in action. See also examples_dir
and get_examples
.
Copied from Trixi.jl.
DispersiveShallowWater.examples_dir
— Methodexamples_dir()
Return the directory where the example files provided with DispersiveShallowWater.jl are located. If DispersiveShallowWater.jl 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())
DispersiveShallowWater.get_examples
— Methodget_examples()
Return a list of all examples that are provided by DispersiveShallowWater.jl. See also examples_dir
and default_example
.
Copied from Trixi.jl.
DispersiveShallowWater.@autoinfiltrate
— Macro@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
Utilities from TrixiBase.jl
Be aware that only trixi_include
is being exported from DispersiveShallowWater.jl. To access the other TrixiBase.jl functions, you need to either:
- Use the fully qualified name:
DispersiveShallowWater.timer()
,DispersiveShallowWater.@trixi_timeit
, etc. - Import TrixiBase.jl explicitly:
using TrixiBase
orimport TrixiBase
TrixiBase.disable_debug_timings
— Methoddisable_debug_timings()
Disable all @trixi_timeit
timings. The timings should be optimized away, allowing for truly zero-overhead. Enable timings again with enable_debug_timings
.
See also enable_debug_timings
, @trixi_timeit
.
TrixiBase.enable_debug_timings
— Methodenable_debug_timings()
Enable all @trixi_timeit
timings (default behavior).
See also disable_debug_timings
, @trixi_timeit
.
TrixiBase.timer
— Methodtimer()
Main timer for global timing, e.g., to be used with @trixi_timeit
.
TrixiBase.trixi_include
— Methodtrixi_include([mapexpr::Function=identity,] [mod::Module=Main,] elixir::AbstractString; kwargs...)
include
the file elixir
and evaluate its content in the global scope of module mod
. You can override specific assignments in elixir
by supplying keyword arguments. Its basic purpose is to make it easier to modify some parameters while running simulations from the REPL. Additionally, this is used in tests to reduce the computational burden for CI while still providing examples with sensible default values for users.
Before replacing assignments in elixir
, the keyword argument maxiters
is inserted into calls to solve
with it's default value used in the SciML ecosystem for ODEs, see the "Miscellaneous" section of the documentation.
The optional first argument mapexpr
can be used to transform the included code before it is evaluated: for each parsed expression expr
in elixir
, the include
function actually evaluates mapexpr(expr)
. If it is omitted, mapexpr
defaults to identity
.
Examples
julia> using TrixiBase, Trixi
julia> redirect_stdout(devnull) do
trixi_include(@__MODULE__, joinpath(examples_dir(), "tree_1d_dgsem", "elixir_advection_extended.jl"),
tspan=(0.0, 0.1))
sol.t[end]
end
[ Info: You just called `trixi_include`. Julia may now compile the code, please be patient.
0.1
TrixiBase.trixi_include_changeprecision
— Methodtrixi_include_changeprecision(T, [mod::Module=Main,] elixir::AbstractString; kwargs...)
include
the elixir elixir
and evaluate its content in the global scope of module mod
. You can override specific assignments in elixir
by supplying keyword arguments, similar to trixi_include
.
The only difference to trixi_include
is that the precision of floating-point numbers in the included elixir is changed to T
. More precisely, the package ChangePrecision.jl is used to convert all Float64
literals, operations like /
that produce Float64
results, and functions like ones
that return Float64
arrays by default, to the desired type T
. See the documentation of ChangePrecision.jl for more details.
The purpose of this function is to conveniently run a full simulation with Float32
, which is orders of magnitude faster on most GPUs than Float64
, by just including the elixir with trixi_include_changeprecision(Float32, elixir)
. Many constructors in the Trixi.jl framework are written in a way that changing all floating-point arguments to Float32
will change the element type to Float32
as well. In TrixiParticles.jl, including an elixir with this macro should be sufficient to run the full simulation with single precision.
TrixiBase.@trixi_timeit
— Macro@trixi_timeit timer() "some label" expression
Basically the same as a special case of @timeit_debug
from TimerOutputs.jl, but without try ... finally ... end
block. Thus, it's not exception-safe, but it also avoids some related performance problems. Since we do not use exception handling in Trixi.jl, that's not really an issue.
All @trixi_timeit
timings can be disabled with disable_debug_timings
. The timings should then be optimized away, allowing for truly zero-overhead.
See also disable_debug_timings
, enable_debug_timings
.