JumpSystem
System Constructors
ModelingToolkit.JumpSystem — Typestruct JumpSystem{U<:RecursiveArrayTools.ArrayPartition} <: ModelingToolkit.AbstractSystemA system of jump processes.
Fields
eqsThe jumps of the system. Allowable types are
ConstantRateJump,VariableRateJump,MassActionJump.
ivThe independent variable, usually time.
statesThe dependent variables, representing the state of the system.
psThe parameters of the system.
observednameThe name of the system.
systemsThe 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.
Example
using ModelingToolkit
@parameters β γ t
@variables S I R
rate₁ = β*S*I
affect₁ = [S ~ S - 1, I ~ I + 1]
rate₂ = γ*I
affect₂ = [I ~ I - 1, R ~ R + 1]
j₁ = ConstantRateJump(rate₁,affect₁)
j₂ = ConstantRateJump(rate₂,affect₂)
j₃ = MassActionJump(2*β+γ, [R => 1], [S => 1, R => -1])
js = JumpSystem([j₁,j₂,j₃], t, [S,I,R], [β,γ])Composition and Accessor Functions
sys.eqsorequations(sys): The equations that define the jump system.sys.statesorstates(sys): The set of states in the jump system.sys.parametersorparameters(sys): The parameters of the jump system.sys.ivorindependent_variable(sys): The independent variable of the jump system.
Problem Constructors
SciMLBase.DiscreteProblem — Typefunction DiffEqBase.DiscreteProblem(sys::JumpSystem, u0map, tspan,
parammap=DiffEqBase.NullParameters; kwargs...)Generates a blank DiscreteProblem for a pure jump JumpSystem to utilize as its prob.prob. This is used in the case where there are no ODEs and no SDEs associated with the system.
Continuing the example from the JumpSystem definition:
using DiffEqBase, DiffEqJump
u₀map = [S => 999, I => 1, R => 0]
parammap = [β => .1/1000, γ => .01]
tspan = (0.0, 250.0)
dprob = DiscreteProblem(js, u₀map, tspan, parammap)DiffEqJump.JumpProblem — Typefunction DiffEqBase.JumpProblem(js::JumpSystem, prob, aggregator; kwargs...)Generates a JumpProblem from a JumpSystem.
Continuing the example from the DiscreteProblem definition:
jprob = JumpProblem(js, dprob, Direct())
sol = solve(jprob, SSAStepper())