ODESystem
System Constructors
ModelingToolkit.ODESystem — Typestruct ODESystem <: ModelingToolkit.AbstractODESystemA system of ordinary differential equations.
Fields
tag: 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.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.tspan: Time span.var_to_name: Array variables.ctrls: Control parameters (some subset ofps).observed: Observed states.tgrad: Time-derivative matrix. Note: this field will not be defined untilcalculate_tgradis called on the system.
jac: Jacobian matrix. Note: this field will not be defined untilcalculate_jacobianis called on the system.
ctrl_jac: Control Jacobian matrix. Note: this field will not be defined untilcalculate_control_jacobianis called on the system.
Wfact:Wfactmatrix. Note: this field will not be defined untilgenerate_factorized_Wis called on the system.
Wfact_t:Wfact_tmatrix. Note: this field will not be defined untilgenerate_factorized_Wis 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 inODEProblem.
torn_matching: torn_matching: Tearing result specifying how to solve the system.
connector_type: connector_type: type of the system
preface: preface: inject assignment statements before the evaluation of the RHS function.
continuous_events: continuous_events: AVector{SymbolicContinuousCallback}that model events. The integrator will use root finding to guarantee that it steps at each zero crossing.
discrete_events: discrete_events: AVector{SymbolicDiscreteCallback}that models events. Symbolic analog toSciMLBase.DiscreteCallbackthat executes an affect when a given condition is true at the end of an integration step.
metadata: metadata: metadata for the system, to be used by downstream packages.
gui_metadata: gui_metadata: metadata for MTK GUI.
tearing_state: tearing_state: cache for intermediate tearing state
substitutions: substitutions: substitutions generated by tearing.
complete: complete: if a modelsysis complete, thensys.xno longer performs namespacing.
discrete_subsystems: discrete_subsystems: a list of discrete subsystems.
unknown_states: unknown_states: a list of actual states needed to be solved by solvers. Only used for ODAEProblem.
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))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.get_u0_p(sys, u0map, parammap)Numeric arrays for the initial condition and parameters givenvar => valuemaps.
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; kwargs...) -> ODESystemPerform the Pantelides algorithm to transform a higher index DAE to an index 1 DAE. kwargs are forwarded to pantelides!. End users are encouraged to call structural_simplify instead, which calls this function internally.
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 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(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.
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_sparsityStandard Problem Constructors
SciMLBase.ODEFunction — MethodDiffEqBase.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 — MethodDiffEqBase.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.SteadyStateProblem — MethodSciMLBase.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 for ODAEProblem(sys::ModelingToolkit.AbstractODESystem, args...). Check Documenter's build log for details.
Expression Constructors
ModelingToolkit.ODEFunctionExpr — TypeODEFunctionExpr{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.DAEFunctionExpr — TypeDAEFunctionExpr{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.SteadyStateProblemExpr — TypeSciMLBase.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.