ODESystem
System Constructors
ModelingToolkit.ODESystem — Typestruct ODESystem <: ModelingToolkit.AbstractODESystemA system of ordinary differential equations.
Fields
eqsThe ODEs defining the system.
ivIndependent variable.
statesDependent (state) variables.
psParameter variables.
observedtgradTime-derivative matrix. Note: this field will not be defined until
calculate_tgradis called on the system.
jacJacobian matrix. Note: this field will not be defined until
calculate_jacobianis called on the system.
WfactWfactmatrix. Note: this field will not be defined untilgenerate_factorized_Wis called on the system.
Wfact_tWfact_tmatrix. Note: this field will not be defined untilgenerate_factorized_Wis called on the system.
nameName: the name of the system
systemssystems: The internal systems
default_u0default_u0: The default initial conditions to use when initial conditions are not supplied in
ODEProblem.
default_pdefault_p: The default parameters to use when parameters are not supplied in
ODEProblem.
structurestructure: structural information of the system
Example
using ModelingToolkit
@parameters t σ ρ β
@variables x(t) y(t) z(t)
D = Differential(t)
eqs = [D(x) ~ σ*(y-x),
D(y) ~ x*(ρ-z)-y,
D(z) ~ x*y - β*z]
de = ODESystem(eqs,t,[x,y,z],[σ,ρ,β])Composition and Accessor Functions
sys.eqsorequations(sys): The equations that define the ODE.sys.statesorstates(sys): The set of states in the ODE.sys.parametersorparameters(sys): The parameters of the ODE.sys.ivorindependent_variable(sys): The independent variable of the ODE.
Transformations
ModelingToolkit.ode_order_lowering — Functionode_order_lowering(sys::ODESystem) -> ODESystem
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_transform — Functionliouville_transform(sys::Any) -> 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 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
Applicable Calculation and Generation Functions
calculate_jacobian
calculate_tgrad
calculate_factorized_W
generate_jacobian
generate_tgrad
generate_factorized_W
jacobian_sparsityProblem Constructors
SciMLBase.ODEFunction — Typefunction 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.ODEProblem — Typefunction 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 iipGenerates an ODEProblem from an ODESystem and allows for automatically symbolically calculating numerical enhancements.