API Index

API Documentation

EarthSciMLBase.CoupledSystemType

A system for composing together other systems using the couple function.

  • systems: Model components to be composed together

  • domaininfo: Initial and boundary conditions and other domain information

  • pdefunctions: A vector of functions where each function takes as an argument the resulting PDESystem after DomainInfo is added to this system, and returns a transformed PDESystem.

  • ops: A vector of operators to run during simulations.
  • callbacks: A vector of callbacks to run during simulations.

  • init_callbacks: Objects x with an init_callback(x, Simulator)::DECallback method.

Things that can be added to a CoupledSystem: * ModelingToolkit.ODESystems. If the ODESystem has a field in the metadata called :coupletype (e.g. ModelingToolkit.get_metadata(sys)[:coupletype] returns a struct type with a single field called sys) then that type will be used to check for methods of EarthSciMLBase.couple that use that type. * Operators * DomainInfos * Callbacks * Types X that implement a EarthSciMLBase.get_callback(::X, ::Simulator)::DECallback method * Other CoupledSystems * Types X that implement a EarthSciMLBase.couple(::X, ::CoupledSystem) or EarthSciMLBase.couple(::CoupledSystem, ::X) method. * Tuples or AbstractVectors of any of the things above.

source
EarthSciMLBase.DomainInfoType

Domain information for a ModelingToolkit.jl PDESystem. It can be used with the + operator to add initial and boundary conditions and coordinate transforms to a ModelingToolkit.jl ODESystem or Catalyst.jl ReactionSystem.

NOTE: The independent variable (usually time) must be first in the list of initial and boundary conditions.

  • partial_derivative_funcs: Function that returns spatial derivatives of the partially-independent variables, optionally performing a coordinate transformation first.

    Current function options in this package are:

    • partialderivatives_δxyδlonlat: Returns partial derivatives after transforming any variables named lat and lon

    from degrees to cartesian meters, assuming a spherical Earth.

    Other packages may implement additional functions. They are encouraged to use function names starting with partialderivatives_.

  • icbc: The sets of initial and/or boundary conditions.
source
EarthSciMLBase.SimulatorType

Specify a simulator for large-scale model runs. Δs represent the grid cell spacing in each dimension; for example Δs = [0.1, 0.1, 1] would represent a grid with 0.1 spacing in the first two dimensions and 1 in the third, in whatever units the grid is natively in. The grid spacings should be given in the same order as the partial independent variables are in the provided DomainInfo.

  • sys::CoupledSystem: The system to be integrated

  • sys_mtk::ModelingToolkit.ODESystem: The ModelingToolkit version of the system

  • domaininfo::DomainInfo: Information about the spatiotemporal simulation domain

  • p::Vector: The system parameter values

  • u_init::Vector: The initial values of the system state variables

  • pvidx::Vector{Int64}: The indexes of the partial independent variables in the system parameter value vector

  • grid::Any: The discretized values of the partial independent variables

  • Δs::Tuple{T, T, T} where T: The spacings for each grid dimension

  • obs_fs::Any: Functions to get the current values of the observed variables with input arguments of time and the partial independent variables

  • obs_fs_idx::Dict{Symbolics.Num, Int64}: Indexes for the obs_fs functions

  • tf_fs::Any: Functions to get the current values of the coordinate transforms with input arguments of time and the partial independent variables

source
EarthSciMLBase.SimulatorIMEXType

A simulator strategy based on implicit-explicit (IMEX) time integration. See here for additional information.

alg is the ODE solver algorithm to use, which should be chosen from the Split ODE Solver algorithms listed here. In most cases, it is recommended to use a Jacobian-free Newton-Krylov linear solution method as described here.

Caution

This strategy does not currently work. Do not use. See here for more details.

  • alg
source
EarthSciMLBase.SimulatorStrangType

A simulator strategy based on Strang splitting. Choose either SimulatorStrangThreads or SimulatorStrangSerial to run the simulation.

Warning

SimulatorStrang strategies will still work if no operator is included, but any callbacks included in the system are executed together with the operators, so if there are no operators in the system, the callbacks will not be executed.

source
EarthSciMLBase.SimulatorStrangSerialType
# Specify the stiff and nonstiff ODE solver algorithm.
# `timestep` is the length of time for each splitting step.
SimulatorStrangSerial(stiffalg, nonstiffalg, timestep)

Perform a simulation using Strang splitting, where the MTK system is assumed to be stiff and the operators are assumed to be non-stiff. The solution will be calculated in serial.

  • stiffalg: Stiff solver algorithm to use (see https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/)

  • nonstiffalg: Non-stiff solver algorithm to use (see https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/)

  • timestep: Length of each splitting time step

source
EarthSciMLBase.SimulatorStrangThreadsType
# Specify the number of threads and the stiff and nonstiff ODE solver algorithm.
# `timestep` is the length of time for each splitting step.
SimulatorStrangThreads(threads, stiffalg, nonstiffalg, timestep)
# Use the default number of threads.
SimulatorStrangThreads(stiffalg, nonstiffalg, timestep)

Perform a simulation using Strang splitting, where the MTK system is assumed to be stiff and the operators are assumed to be non-stiff. The solution of the stiff ODE system is parallelized across grid cells using the specified number of threads.

  • threads: Number of threads to use

  • stiffalg: Stiff solver algorithm to use (see https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/)

  • nonstiffalg: Non-stiff solver algorithm to use (see https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/)

  • timestep: Length of each splitting time step

source
EarthSciMLBase.SimulatorStrategyType

SimulatorStrategy is an abstract type that defines the strategy for running a simulation. Each SimulatorStrategy should implement a method of:

run!(st::SimulatorStrategy, s::Simulator{T}, u0) where T

where u0 is the initial conditions for the system state.

source
EarthSciMLBase.constBCType

Construct constant boundary conditions equal to the value specified by val.

  • val: The value of the constant boundary conditions.

  • partialdomains: The partial domains, e.g. [x ∈ Interval(x_min, x_max), y ∈ Interval(y_min, y_max)].

source
EarthSciMLBase.constICType

Construct constant initial conditions equal to the value specified by val.

  • val: The value of the constant initial conditions.

  • indepdomain: The independent domain, e.g. t ∈ Interval(t_min, t_max).

source
EarthSciMLBase.periodicBCType

Construct periodic boundary conditions for the given partialdomains. Periodic boundary conditions are defined as when the value at one side of the domain is set equal to the value at the other side, so that the domain "wraps around" from one side to the other.

  • partialdomains: The partial domains, e.g. [x ∈ Interval(x_min, x_max), y ∈ Interval(y_min, y_max)].
source
EarthSciMLBase.zerogradBCType

Construct zero-gradient boundary conditions for the given partialdomains.

  • partialdomains: The partial domains, e.g. [x ∈ Interval(x_min, x_max), y ∈ Interval(y_min, y_max)].
source
EarthSciMLBase.ConstantWindMethod
ConstantWind(t, vals; name)

Construct a constant wind velocity model component with the given wind speed(s), which should include units. For example, ConstantWind(t, 1u"m/s", 2u"m/s").

source
EarthSciMLBase.MeanWindMethod
MeanWind(t, domain)

A model component that represents the mean wind velocity, where pvars is the partial dependent variables for the domain.

source
EarthSciMLBase.add_dimsMethod
add_dims(expression, vars, dims)
add_dims(equation, vars, dims)

Add the given dimensions to each variable in vars in the given expression or equation. Each variable in vars must be unidimensional, i.e. defined like @variables u(t) rather than @variables u(..).

Example:

using EarthSciMLBase, ModelingToolkit

@parameters x y k t
@variables u(t) q(t)
exp = 2u + 3k*q + 1
EarthSciMLBase.add_dims(exp, [u, q], [x, y, t])

# output
1 + 2u(x, y, t) + 3k*q(x, y, t)
source
EarthSciMLBase.add_scopeMethod
add_scope(sys, v, iv)

Add a system scope to a variable name, for example so that x in system sys1 becomes sys1₊x. iv is the independent variable.

source
EarthSciMLBase.coupleMethod
couple(systems...) -> CoupledSystem

Couple multiple ModelingToolkit systems together.

The systems that are arguments to this system can be of type ModelingToolkit.AbstractSystem, CoupledSystem, DomainInfo, or any type T that has a method couple(::CoupledSystem, ::T)::CoupledSystem or a method couple(::T, ::CoupledSystem)::CoupledSystem defined for it.

source
EarthSciMLBase.couple2Method
couple2()

Perform bi-directional coupling for two equation systems.

To specify couplings for system pairs, create methods for this function with the signature:

EarthSciMLBase.couple2(a::ACoupler, b::BCoupler)::ConnectorSystem

where ACoupler and BCoupler are :coupletypes defined like this:

struct ACoupler sys end
@named asys = ODESystem([], t, metadata=Dict(:coupletype=>ACoupler))
source
EarthSciMLBase.dimsMethod
dims(
    icbc::EarthSciMLBase.ICcomponent
) -> Vector{Symbolics.Num}

Returns the dimensions of the independent and partial domains associated with these initial or boundary conditions.

source
EarthSciMLBase.domainsMethod
domains(icbc::EarthSciMLBase.ICcomponent) -> Vector

Returns the domains associated with these initial or boundary conditions.

source
EarthSciMLBase.get_dvMethod

Return the dependent variable, which is the first argument of the term, unless the term is a time derivative, in which case the dependent variable is the argument of the time derivative.

source
EarthSciMLBase.get_needed_varsMethod
get_needed_vars(sys)

Return the indexes of the system variables that the state variables of the final simplified system depend on. This should be done before running structural_simplify on the system.

source
EarthSciMLBase.gridMethod
grid(d, Δs)

Return the ranges representing the discretization of the partial independent variables for this domain, based on the discretization intervals given in Δs

source
EarthSciMLBase.icbcMethod
icbc(di, states)

Return a vector of equations that define the initial and boundary conditions for the given state variables.

source
EarthSciMLBase.init_callbackMethod

Types that implement an:

init_callback(x, Simulator)::DECallback

method can also be coupled into a CoupledSystem. The init_callback function will be run before the simulator is run to get the callback.

source
EarthSciMLBase.ivarMethod
ivar(di::DomainInfo) -> Any

Return the independent variable associated with these initial and boundary conditions.

source
EarthSciMLBase.mtk_opMethod

Return a SciMLOperator to apply the MTK system to each column of s.u after reshaping to a matrix.

source
EarthSciMLBase.observed_expressionMethod
observed_expression(eqs, x)

Return an expression for the observed value of a variable x after substituting in the constants observed values of other variables. extra_eqs is a list of additional equations to use in the substitution.

source
EarthSciMLBase.observed_functionMethod
observed_function(eqs, x, coords)

Return a function to for the observed value of a variable x based on the input arguments in coords. extra_eqs is a list of additional equations to use to determine the value of x.

source
EarthSciMLBase.operator_composeFunction
operator_compose(a, b)
operator_compose(a, b, translate)

Compose to systems of equations together by adding the right-hand side terms together of equations that have matching left-hand sides. The left hand sides of two equations will be considered matching if:

  1. They are both time derivatives of the same variable.
  2. The first one is a time derivative of a variable and the second one is the variable itself.
  3. There is an entry in the optional translate dictionary that maps the dependent variable in the first system to the dependent variable in the second system, e.g. Dict(sys1.sys.x => sys2.sys.y).
  4. There is an entry in the optional translate dictionary that maps the dependent variable in the first system to the dependent variable in the second system, with a conversion factor, e.g. Dict(sys1.sys.x => sys2.sys.y => 6).
source
EarthSciMLBase.param_to_varMethod

Replace the parameter p in the system sys with a new variable that has the same name, units, and description as p.

param_to_var(sys, ps)

This can be useful to replace a parameter that does not change in time in a model component with one specified by another system that does change in time (or space). For example, the code below specifies a first-order loss equation, and then changes the temperature (which determines the loss rate) with a temperature value that varies in time. ```

source
EarthSciMLBase.partialderivatives_δxyδlonlatMethod
partialderivatives_δxyδlonlat(pvars; default_lat)

Return partial derivative operator transform factors corresponding for the given partial-independent variables after converting variables named lon and lat from degrees to x and y meters, assuming they represent longitude and latitude on a spherical Earth.

source
EarthSciMLBase.prune_observedMethod
prune_observed(sys)

Remove equations from an ODESystem where the variable in the LHS is not present in any of the equations for the state variables. This can be used to remove computationally intensive equations that are not used in the final model.

source
EarthSciMLBase.pvarsMethod
pvars(di::DomainInfo) -> Any

Return the partial independent variables associated with these initial and boundary conditions.

source
EarthSciMLBase.run!Function
run!(s::Simulator, st::SimulatorIMEX; ...)
run!(s::Simulator, st::SimulatorIMEX, u; kwargs...)

Run the simulation. kwargs are passed to the ODEProblem and ODE solver constructors.

source
EarthSciMLBase.run!Method
run!(
    s::Simulator{T},
    st::EarthSciMLBase.SimulatorStrang;
    ...
) -> Any
run!(
    s::Simulator{T},
    st::EarthSciMLBase.SimulatorStrang,
    u;
    kwargs...
) -> Any

Run the simualation. kwargs are passed to the ODEProblem and ODE solver constructors.

source
EarthSciMLBase.steplengthMethod
steplength(timesteps)

Return the time step length common to all of the given timesteps. Throw an error if not all timesteps are the same length.

source
EarthSciMLBase.utypeMethod
utype(_)

Return the data type of the state variables for this domain, based on the data types of the boundary conditions domain intervals.

source