ODESystem

System Constructors

ModelingToolkit.ODESystemType
struct ODESystem <: ModelingToolkit.AbstractODESystem

A system of ordinary differential equations.

Fields

  • tag: A tag for the system. If two systems have the same tag, then they are structurally identical.
  • eqs: The ODEs defining the system.

  • iv: Independent variable.

  • unknowns: Dependent (unknown) variables. Must not contain the independent variable.

    N.B.: If torn_matching !== nothing, this includes all variables. Actual ODE unknowns are determined by the SelectedState() entries in torn_matching.

  • ps: Parameter variables. Must not contain the independent variable.

  • tspan: Time span.

  • var_to_name: Array variables.

  • ctrls: Control parameters (some subset of ps).

  • observed: Observed equations.

  • constraintsystem: System of constraints that must be satisfied by the solution to the system.

  • costs: A set of expressions defining the costs of the system for optimal control.

  • consolidate: Takes the cost vector and returns a scalar for optimization.

  • tgrad: Time-derivative matrix. Note: this field will not be defined until calculate_tgrad is called on the system.

  • jac: Jacobian matrix. Note: this field will not be defined until calculate_jacobian is called on the system.
  • name: The name of the system.
  • description: A description of the system.
  • systems: The internal systems. These are required to have unique names.
  • defaults: The default values to use when initial conditions and/or parameters are not supplied in ODEProblem.
  • guesses: The guesses to use as the initial conditions for the initialization system.
  • torn_matching: Tearing result specifying how to solve the system.
  • initializesystem: The system for performing the initialization.
  • initialization_eqs: Extra equations to be enforced during the initialization sequence.
  • schedule: The schedule for the code generation process.
  • connector_type: Type of the system.
  • preface: Inject assignment statements before the evaluation of the RHS function.
  • continuous_events: A Vector{SymbolicContinuousCallback} that model events. The integrator will use root finding to guarantee that it steps at each zero crossing.
  • discrete_events: A Vector{SymbolicDiscreteCallback} that models events. Symbolic analog to SciMLBase.DiscreteCallback that executes an affect when a given condition is true at the end of an integration step.
  • parameter_dependencies: Topologically sorted parameter dependency equations, where all symbols are parameters and the LHS is a single parameter.
  • assertions: Mapping of conditions which should be true throughout the solution process to corresponding error messages. These will be added to the equations when calling debug_system.
  • metadata: Metadata for the system, to be used by downstream packages.
  • gui_metadata: Metadata for MTK GUI.
  • is_dde: A boolean indicating if the given ODESystem represents a system of DDEs.
  • tstops: A list of points to provide to the solver as tstops. Uses the same syntax as discrete events.
  • tearing_state: Cache for intermediate tearing state.
  • substitutions: Substitutions generated by tearing.
  • namespacing: If false, then sys.x no longer performs namespacing.
  • complete: If true, denotes the model will not be modified any further.
  • index_cache: Cached data for fast symbolic indexing.
  • discrete_subsystems: A list of discrete subsystems.
  • solved_unknowns: A list of actual unknowns needed to be solved by solvers.
  • split_idxs: A vector of vectors of indices for the split parameters.
  • ignored_connections: The analysis points removed by transformations, representing connections to be ignored. The first element of the tuple analysis points connecting systems and the second are ones connecting variables (for the trivial form of connect).
  • parent: The hierarchical parent system before simplification.

Example

using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D

@parameters σ ρ β
@variables x(t) y(t) z(t)

eqs = [D(x) ~ σ*(y-x),
       D(y) ~ x*(ρ-z)-y,
       D(z) ~ x*y - β*z]

@named de = ODESystem(eqs,t,[x,y,z],[σ,ρ,β],tspan=(0, 1000.0))
source

Composition and Accessor Functions

  • get_eqs(sys) or equations(sys): The equations that define the ODE.
  • get_unknowns(sys) or unknowns(sys): The set of unknowns in the ODE.
  • get_ps(sys) or parameters(sys): The parameters of the ODE.
  • get_iv(sys): The independent variable of the ODE.
  • get_u0_p(sys, u0map, parammap) Numeric arrays for the initial condition and parameters given var => value maps.
  • continuous_events(sys): The set of continuous events in the ODE.
  • discrete_events(sys): The set of discrete events in the ODE.
  • alg_equations(sys): The algebraic equations (i.e. that does not contain a differential) that defines the ODE.
  • get_alg_eqs(sys): The algebraic equations (i.e. that does not contain a differential) that defines the ODE. Only returns equations of the current-level system.
  • diff_equations(sys): The differential equations (i.e. that contain a differential) that defines the ODE.
  • get_diff_eqs(sys): The differential equations (i.e. that contain a differential) that defines the ODE. Only returns equations of the current-level system.
  • has_alg_equations(sys): Returns true if the ODE contains any algebraic equations (i.e. that does not contain a differential).
  • has_alg_eqs(sys): Returns true if the ODE contains any algebraic equations (i.e. that does not contain a differential). Only considers the current-level system.
  • has_diff_equations(sys): Returns true if the ODE contains any differential equations (i.e. that does contain a differential).
  • has_diff_eqs(sys): Returns true if the ODE contains any differential equations (i.e. that does contain a differential). Only considers the current-level system.

Transformations

ModelingToolkit.structural_simplifyFunction
structural_simplify(sys; ...)
structural_simplify(
    sys,
    io;
    additional_passes,
    simplify,
    split,
    allow_symbolic,
    allow_parameter,
    conservative,
    fully_determined,
    kwargs...
)

Structurally simplify algebraic equations in a system and compute the topological sort of the observed equations in sys.

Optional Arguments:

  • optional argument io may take a tuple (inputs, outputs). This will convert all inputs to parameters and allow them to be unconnected, i.e., simplification will allow models where n_unknowns = n_equations - n_inputs.

Optional Keyword Arguments:

  • When simplify=true, the simplify function will be applied during the tearing process.
  • allow_symbolic=false, allow_parameter=true, and conservative=false limit the coefficient types during tearing. In particular, conservative=true limits tearing to only solve for trivial linear systems where the coefficient has the absolute value of $1$.
  • fully_determined=true controls whether or not an error will be thrown if the number of equations don't match the number of inputs, outputs, and equations.
  • sort_eqs=true controls whether equations are sorted lexicographically before simplification or not.
source
ModelingToolkit.ode_order_loweringFunction
ode_order_lowering(sys::ODESystem) -> Any

Takes a Nth order ODESystem and returns a new ODESystem written in first order form by defining new variables which represent the N-1 derivatives.

source
ModelingToolkit.change_independent_variableFunction
change_independent_variable(
    sys::AbstractODESystem, iv, eqs = [];
    add_old_diff = false, simplify = true, fold = false
)

Transform the independent variable (e.g. $t$) of the ODE system sys to a dependent variable iv (e.g. $u(t)$). The transformation is well-defined when the mapping between the new and old independent variables are one-to-one. This is satisfied if one is a strictly increasing function of the other (e.g. $du(t)/dt > 0$ or $du(t)/dt < 0$).

Any extra equations eqs involving the new and old independent variables will be taken into account in the transformation.

Keyword arguments

  • add_old_diff: Whether to add a differential equation for the old independent variable in terms of the new one using the inverse function rule $dt/du = 1/(du/dt)$.
  • simplify: Whether expanded derivative expressions are simplified. This can give a tidier transformation.
  • fold: Whether internal substitutions will evaluate numerical expressions.

Usage before structural simplification

The variable change must take place before structural simplification. In following calls to structural_simplify, consider passing allow_symbolic = true to avoid undesired constraint equations between between dummy variables.

Usage with non-autonomous systems

If sys is non-autonomous (i.e. $t$ appears explicitly in its equations), consider passing an algebraic equation relating the new and old independent variables (e.g. $t = f(u(t))$). Otherwise the transformed system can be underdetermined. If an algebraic relation is not known, consider using add_old_diff instead.

Usage with hierarchical systems

It is recommended that iv is a non-namespaced variable in sys. This means it can belong to the top-level system or be a variable in a subsystem declared with GlobalScope.

Example

Consider a free fall with constant horizontal velocity. Physics naturally describes position as a function of time. By changing the independent variable, it can be reformulated for vertical position as a function of horizontal position:

julia> @variables x(t) y(t);

julia> @named M = ODESystem([D(D(y)) ~ -9.81, D(D(x)) ~ 0.0], t);

julia> M = change_independent_variable(M, x);

julia> M = structural_simplify(M; allow_symbolic = true);

julia> unknowns(M)
3-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 xˍt(x)
 y(x)
 yˍx(x)
source
ModelingToolkit.liouville_transformFunction
liouville_transform(
    sys::ModelingToolkit.AbstractODESystem;
    kwargs...
) -> ODESystem

Generates the Liouville transformed set of ODEs, which is the original ODE system with a new variable trJ appended, corresponding to the -tr(Jacobian). This variable is used for properties like uncertainty propagation from a given initial distribution density.

For example, if $u'=p*u$ and p follows a probability distribution $f(p)$, then the probability density of a future value with a given choice of $p$ is computed by setting the initial trJ = f(p), and the final value of trJ is the probability of $u(t)$.

Example:

using ModelingToolkit, OrdinaryDiffEq

@independent_variables t
@parameters α β γ δ
@variables x(t) y(t)
D = Differential(t)
eqs = [D(x) ~ α*x - β*x*y, D(y) ~ -δ*y + γ*x*y]
@named sys = ODESystem(eqs, t)

sys2 = liouville_transform(sys)
sys2 = complete(sys2)
u0 = [x => 1.0, y => 1.0, sys2.trJ => 1.0]
prob = ODEProblem(sys2, u0, tspan, p)
sol = solve(prob, Tsit5())

Where sol[3,:] is the evolution of trJ over time.

Sources:

Probabilistic Robustness Analysis of F-16 Controller Performance: An Optimal Transport Approach

Abhishek Halder, Kooktae Lee, and Raktim Bhattacharya https://abhishekhalder.bitbucket.io/F16ACC2013Final.pdf

source
Missing docstring.

Missing docstring for alias_elimination. Check Documenter's build log for details.

Analyses

Missing docstring.

Missing docstring for ModelingToolkit.islinear. Check Documenter's build log for details.

Missing docstring.

Missing docstring for ModelingToolkit.isautonomous. Check Documenter's build log for details.

Missing docstring.

Missing docstring for ModelingToolkit.isaffine. Check Documenter's build log for details.

Applicable Calculation and Generation Functions

ModelingToolkit.calculate_jacobianFunction
calculate_jacobian(sys::AbstractSystem)

Calculate the Jacobian matrix of a system.

Returns a matrix of Num instances. The result from the first call will be cached in the system object.

source
ModelingToolkit.calculate_tgradFunction
calculate_tgrad(sys::AbstractTimeDependentSystem)

Calculate the time gradient of a system.

Returns a vector of Num instances. The result from the first call will be cached in the system object.

source
ModelingToolkit.calculate_factorized_WFunction
calculate_factorized_W(sys::AbstractSystem)

Calculate the factorized W-matrix of a system.

Returns a matrix of Num instances. The result from the first call will be cached in the system object.

source
ModelingToolkit.generate_jacobianFunction
generate_jacobian(sys::AbstractSystem, dvs = unknowns(sys), ps = parameters(sys),
                  expression = Val{true}; sparse = false, kwargs...)

Generates a function for the Jacobian matrix of a system. Extra arguments control the arguments to the internal build_function call.

source
ModelingToolkit.generate_tgradFunction
generate_tgrad(sys::AbstractTimeDependentSystem, dvs = unknowns(sys), ps = parameters(sys),
               expression = Val{true}; kwargs...)

Generates a function for the time gradient of a system. Extra arguments control the arguments to the internal build_function call.

source
ModelingToolkit.generate_factorized_WFunction
generate_factorized_W(sys::AbstractSystem, dvs = unknowns(sys), ps = parameters(sys),
                      expression = Val{true}; sparse = false, kwargs...)

Generates a function for the factorized W matrix of a system. Extra arguments control the arguments to the internal build_function call.

source
Missing docstring.

Missing docstring for jacobian_sparsity. Check Documenter's build log for details.

Standard Problem Constructors

SciMLBase.ODEFunctionMethod
DiffEqBase.ODEFunction{iip}(sys::AbstractODESystem, dvs = unknowns(sys),
                            ps = parameters(sys);
                            version = nothing, tgrad = false,
                            jac = false,
                            sparse = false,
                            kwargs...) where {iip}

Create an ODEFunction from the ODESystem. The arguments dvs and ps are used to set the order of the dependent variable and parameter vectors, respectively.

source
SciMLBase.ODEProblemMethod
DiffEqBase.ODEProblem{iip}(sys::AbstractODESystem, u0map, tspan,
                           parammap = DiffEqBase.NullParameters();
                           allow_cost = false,
                           version = nothing, tgrad = false,
                           jac = false,
                           checkbounds = false, sparse = false,
                           simplify = false,
                           linenumbers = true, parallel = SerialForm(),
                           kwargs...) where {iip}

Generates an ODEProblem from an ODESystem and allows for automatically symbolically calculating numerical enhancements.

source
SciMLBase.SteadyStateProblemMethod
SciMLBase.SteadyStateProblem(sys::AbstractODESystem, u0map,
                             parammap = DiffEqBase.NullParameters();
                             version = nothing, tgrad = false,
                             jac = false,
                             checkbounds = false, sparse = false,
                             linenumbers = true, parallel = SerialForm(),
                             kwargs...) where {iip}

Generates an SteadyStateProblem from an ODESystem and allows for automatically symbolically calculating numerical enhancements.

source
Missing docstring.

Missing docstring for DAEProblem(sys::ModelingToolkit.AbstractODESystem, args...). Check Documenter's build log for details.

Expression Constructors

ModelingToolkit.ODEFunctionExprType
ODEFunctionExpr{iip}(sys::AbstractODESystem, dvs = unknowns(sys),
                     ps = parameters(sys);
                     version = nothing, tgrad = false,
                     jac = false,
                     sparse = false,
                     kwargs...) where {iip}

Create a Julia expression for an ODEFunction from the ODESystem. The arguments dvs and ps are used to set the order of the dependent variable and parameter vectors, respectively.

source
ModelingToolkit.DAEFunctionExprType
DAEFunctionExpr{iip}(sys::AbstractODESystem, dvs = unknowns(sys),
                     ps = parameters(sys);
                     version = nothing, tgrad = false,
                     jac = false,
                     sparse = false,
                     kwargs...) where {iip}

Create a Julia expression for an ODEFunction from the ODESystem. The arguments dvs and ps are used to set the order of the dependent variable and parameter vectors, respectively.

source
ModelingToolkit.SteadyStateProblemExprType
SciMLBase.SteadyStateProblemExpr(sys::AbstractODESystem, u0map,
                                 parammap = DiffEqBase.NullParameters();
                                 version = nothing, tgrad = false,
                                 jac = false,
                                 checkbounds = false, sparse = false,
                                 skipzeros = true, fillzeros = true,
                                 linenumbers = true, parallel = SerialForm(),
                                 kwargs...) where {iip}

Generates a Julia expression for building a SteadyStateProblem from an ODESystem and allows for automatically symbolically calculating numerical enhancements.

source