ODESystem

System Constructors

ModelingToolkit.ODESystemType
struct ODESystem <: ModelingToolkit.AbstractODESystem

A system of ordinary differential equations.

Fields

  • eqs

    The ODEs defining the system.

  • iv

    Independent variable.

  • states

    Dependent (state) variables. Must not contain the independent variable.

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

  • ps

    Parameter variables. Must not contain the independent variable.

  • var_to_name

    Array variables.

  • ctrls

    Control parameters (some subset of ps).

  • observed

    Observed states.

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

  • Wfact

    Wfact matrix. Note: this field will not be defined until generate_factorized_W is called on the system.

  • Wfact_t

    Wfact_t matrix. Note: this field will not be defined until generate_factorized_W is called on the system.

  • name

    Name: the name of the system

  • systems

    systems: The internal systems. These are required to have unique names.

  • defaults

    defaults: The default values to use when initial conditions and/or parameters are not supplied in ODEProblem.

  • torn_matching

    torn_matching: Tearing result specifying how to solve the system.

  • connector_type

    connector_type: type of the system

  • connections

    connections: connections in a system

  • preface

    preface: inject assignment statements before the evaluation of the RHS function.

  • continuous_events

    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

    discrete_events: A Vector{SymbolicDiscreteCallback} that models events. Symbolic analog to SciMLBase.DiscreteCallback that exectues an affect when a given condition is true at the end of an integration step.

  • tearing_state

    tearing_state: cache for intermediate tearing state

  • substitutions

    substitutions: substitutions generated by tearing.

Example

using ModelingToolkit

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

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

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

Composition and Accessor Functions

  • get_eqs(sys) or equations(sys): The equations that define the ODE.
  • get_states(sys) or states(sys): The set of states in the ODE.
  • get_ps(sys) or parameters(sys): The parameters of the ODE.
  • get_iv(sys): The independent variable of the ODE.

Transformations

Missing docstring.

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

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.

ModelingToolkit.liouville_transformFunction
liouville_transform(sys::ModelingToolkit.AbstractODESystem)

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 inital trJ = f(p), and the final value of trJ is the probability of $u(t)$.

Example:

using ModelingToolkit, OrdinaryDiffEq, Test

@parameters t α β γ δ
@variables x(t) y(t)
D = Differential(t)

eqs = [D(x) ~ α*x - β*x*y,
       D(y) ~ -δ*y + γ*x*y]

sys = ODESystem(eqs)
sys2 = liouville_transform(sys)
@variables trJ

u0 = [x => 1.0,
      y => 1.0,
      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

Missing docstring.

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

Missing docstring.

Missing docstring for tearing. 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

calculate_jacobian
calculate_tgrad
calculate_factorized_W
generate_jacobian
generate_tgrad
generate_factorized_W
jacobian_sparsity

Standard Problem Constructors

SciMLBase.ODEFunctionMethod
function DiffEqBase.ODEFunction{iip}(sys::AbstractODESystem, dvs = states(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.

SciMLBase.ODEProblemMethod
function DiffEqBase.ODEProblem{iip}(sys::AbstractODESystem,u0map,tspan,
                                    parammap=DiffEqBase.NullParameters();
                                    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.

SciMLBase.SteadyStateProblemMethod
function 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.

Torn Problem Constructors

Missing docstring.

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

Expression Constructors

ModelingToolkit.ODEFunctionExprType
function ODEFunctionExpr{iip}(sys::AbstractODESystem, dvs = states(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.

ModelingToolkit.DAEFunctionExprType
function DAEFunctionExpr{iip}(sys::AbstractODESystem, dvs = states(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.

ModelingToolkit.SteadyStateProblemExprType
function 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.