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

  • 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.
  • 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: A mapping from dependent parameters to expressions describing how they are calculated from other parameters.
  • metadata: Metadata for the system, to be used by downstream packages.
  • gui_metadata: Metadata for MTK GUI.
  • tearing_state: Cache for intermediate tearing state.
  • substitutions: Substitutions generated by tearing.
  • complete: If a model sys is complete, then sys.x no longer performs namespacing.
  • 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.
  • parent: The hierarchical parent system before simplification.

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],[σ,ρ,β],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; simplify, split, kwargs...)

Structurally simplify algebraic equations in a system and compute the topological sort of the observed equations. When simplify=true, the simplify function will be applied during the tearing process. It also takes kwargs allow_symbolic=false and allow_parameter=true which limits the coefficient types during tearing.

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

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.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 initial 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(complete(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 = full_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 = full_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 = full_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();
                           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

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