JumpSystem

System Constructors

ModelingToolkit.JumpSystemType
struct JumpSystem{U<:RecursiveArrayTools.ArrayPartition} <: ModelingToolkit.AbstractSystem

A system of jump processes.

Fields

  • eqs

    The jumps of the system. Allowable types are ConstantRateJump, VariableRateJump, MassActionJump.

  • iv

    The independent variable, usually time.

  • states

    The dependent variables, representing the state of the system.

  • ps

    The parameters of the system.

  • pins

  • observed

  • name

    The name of the system.

  • systems

    The internal systems.

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], [β,γ])
source

Composition and Accessor Functions

  • sys.eqs or equations(sys): The equations that define the jump system.
  • sys.states or states(sys): The set of states in the jump system.
  • sys.parameters or parameters(sys): The parameters of the jump system.
  • sys.iv or independent_variable(sys): The independent variable of the jump system.

Problem Constructors

DiffEqBase.DiscreteProblemType
function 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)
source
DiffEqJump.JumpProblemType
function 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())
source