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. Must not contain the independent variable.
psParameter variables. Must not contain the independent variable.
var_to_nameArray variables.
ctrlsControl parameters (some subset of
ps).observedObserved states.
tgradTime-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.
ctrl_jacControl Jacobian matrix. Note: this field will not be defined until
calculate_control_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. These are required to have unique names.
defaultsdefaults: The default values to use when initial conditions and/or parameters are not supplied in
ODEProblem.
structurestructure: structural information of the system
connection_typeconnection_type: type of the system
prefacepreface: injuect assignment statements before the evaluation of the RHS function.
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]
de = ODESystem(eqs,t,[x,y,z],[σ,ρ,β])Composition and Accessor Functions
get_eqs(sys)orequations(sys): The equations that define the ODE.get_states(sys)orstates(sys): The set of states in the ODE.get_ps(sys)orparameters(sys): The parameters of the ODE.get_iv(sys): The independent variable of the ODE.
Transformations
Missing docstring for structural_simplify. Check Documenter's build log for details.
ModelingToolkit.ode_order_lowering — Functionode_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.StructuralTransformations.dae_index_lowering — Functiondae_index_lowering(sys::ODESystem) -> ODESystemPerform the Pantelides algorithm to transform a higher index DAE to an index 1 DAE.
ModelingToolkit.liouville_transform — Functionliouville_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 for alias_elimination. Check Documenter's build log for details.
Missing docstring for tearing. Check Documenter's build log for details.
Analyses
Missing docstring for ModelingToolkit.islinear. Check Documenter's build log for details.
Missing docstring for ModelingToolkit.isautonomous. 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_sparsityStandard Problem 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.
Missing docstring for SteadyStateFunction. Check Documenter's build log for details.
SciMLBase.SteadyStateProblem — Typefunction DiffEqBase.SteadyStateProblem(sys::AbstractODESystem,u0map,
parammap=DiffEqBase.NullParameters();
version = nothing, tgrad=false,
jac = false,
checkbounds = false, sparse = false,
linenumbers = true, parallel=SerialForm(),
kwargs...) where iipGenerates an SteadyStateProblem from an ODESystem and allows for automatically symbolically calculating numerical enhancements.
Torn Problem Constructors
Missing docstring for ODAEProblem. Check Documenter's build log for details.