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

Model specific equations

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 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
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.

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.

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=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.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.energy_totalMethod
energy_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.

source
DispersiveShallowWater.entropyMethod
entropy(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.

source
DispersiveShallowWater.nondim2primMethod
nondim2prim(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.

Parameter constraints

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.

source
DispersiveShallowWater.prim2nondimMethod
prim2nondim(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.

Parameter constraints

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.

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.initial_condition_manufactured_reflectingMethod
initial_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
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

Linear dispersion relations

DispersiveShallowWater.EulerEquations1DType
EulerEquations1D(; gravity, eta0 = 0.0)

A struct representing the 1D Euler equations with a given gravitational acceleration gravity and a still-water surface eta0.

Note

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.

source
DispersiveShallowWater.LinearDispersionRelationType
LinearDispersionRelation(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.

source
DispersiveShallowWater.wave_speedMethod
wave_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.

source

Mesh

Boundary conditions

Solver

DispersiveShallowWater.SolverMethod
Solver(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.

Periodic operators only

This constructor creates periodic derivative operators that are only compatible with periodic boundary conditions. For non-periodic boundary conditions, use the Solver(D1, D2, D3) constructor with appropriate non-periodic operators (or nothing).

source
DispersiveShallowWater.SolverMethod
Solver(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.

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)),
                                    tmp_partitioned = allocate_coefficients(mesh, equations, solver)))

Construct a semidiscretization of a PDE.

source
DispersiveShallowWater.check_solverMethod
check_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.

source
DispersiveShallowWater.jacobianMethod
DispersiveShallowWater.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.

Warning

This functionality may not be implemented for all equations.

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

Experimental data

DispersiveShallowWater.data_dingemansMethod
data_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).

source

Utilities from DispersiveShallowWater.jl

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.data_dirMethod
data_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())
source
DispersiveShallowWater.examples_dirMethod
examples_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())
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

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 or import TrixiBase
TrixiBase.trixi_includeMethod
trixi_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
source
TrixiBase.trixi_include_changeprecisionMethod
trixi_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.

source
TrixiBase.@trixi_timeitMacro
@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.

source