JumpProcesses.jl API
Core Types
JumpProcesses.JumpProblem — Typemutable struct JumpProblem{iip, P, A, C, J<:Union{Nothing, JumpProcesses.AbstractJumpAggregator}, J1, J2, J3, J4, R, K} <: SciMLBase.AbstractJumpProblem{P, J<:Union{Nothing, JumpProcesses.AbstractJumpAggregator}}Defines a collection of jump processes to associate with another problem type.
Constructors
JumpProblems can be constructed by first building another problem type to which the jumps will be associated. For example, to simulate a collection of jump processes for which the transition rates are constant between jumps (called ConstantRateJumps or MassActionJumps), we must first construct a DiscreteProblem
prob = DiscreteProblem(u0, p, tspan)where u0 is the initial condition, p the parameters and tspan the time span. If we wanted to have the jumps coupled with a system of ODEs, or have transition rates with explicit time dependence, we would use an ODEProblem instead that defines the ODE portion of the dynamics.
Given prob we define the jumps via
JumpProblem(prob, aggregator::AbstractAggregatorAlgorithm, jumps::JumpSet ; kwargs...)JumpProblem(prob, aggregator::AbstractAggregatorAlgorithm, jumps...; kwargs...)
Here aggregator specifies the underlying algorithm for calculating next jump times and types, for example Direct. The collection of different AbstractJump types can then be passed within a single JumpSet or as subsequent sequential arguments.
Fields
prob: The type of problem to couple the jumps to. For a pure jump process useDiscreteProblem, to couple to ODEs,ODEProblem, etc.aggregator: The aggregator algorithm that determines the next jump times and types forConstantRateJumps andMassActionJumps. Examples includeDirect.discrete_jump_aggregation: The underlying state data associated with the chosen aggregator.jump_callback:CallBackSetwith the underlyingConstantRateandVariableRatejumps.constant_jumps: TheConstantRateJumps.variable_jumps: TheVariableRateJumps.regular_jump: TheRegularJumps.massaction_jump: TheMassActionJumps.rng: The random number generator to use.kwargs: kwargs to pass on to solve call.
Keyword Arguments
rng, the random number generator to use. Defaults to Julia's built-in generator.save_positions=(true,true)when including variable rates and(false,true)for constant rates, specifies whether to save the system's state (before, after) the jump occurs.spatial_system, for spatial problems the underlying spatial structure.hopping_constants, for spatial problems the spatial transition rate coefficients.use_vrj_bounds = true, set to false to disable handling boundedVariableRateJumps with a supporting aggregator (such asCoevolve). They will then be handled via the continuous integration interface, and treated like generalVariableRateJumps.vr_aggregator, indicates the aggregator to use for sampling variable rate jumps. Current default isVR_FRM.
Please see the tutorial page in the DifferentialEquations.jl docs for usage examples and commonly asked questions.
JumpProcesses.SSAStepper — Typestruct SSAStepper <: SciMLBase.AbstractDEAlgorithmHighly efficient integrator for pure jump problems that involve only ConstantRateJumps, MassActionJumps, and/or VariableRateJumps with rate bounds.
Notes
- Only works with
JumpProblems defined fromDiscreteProblems. - Only works with collections of
ConstantRateJumps,MassActionJumps, andVariableRateJumps with rate bounds. - Only supports
DiscreteCallbacks for events, which are checked after every step taken bySSAStepper. - Only supports a limited subset of the output controls from the common solver interface, specifically
save_start,save_end, andsaveat. - As when using jumps with ODEs and SDEs, saving controls for whether to save each time a jump occurs are via the
save_positionskeyword argument toJumpProblem. Note that when choosingSSAStepperas the timestepper,save_positions = (true,true),(true,false), or(false,true)are all equivalent.SSAStepperwill save only the post-jump state in the solution object in each of these cases. This is because solution objects generated viaSSAStepperuse piecewise-constant interpolation, and can therefore exactly reconstruct the sampled jump process path with knowing just the post-jump state. That is,sol(t)for any0 <= t <= tstopwill give the exact value of the sampled solution path attprovided at least one component ofsave_positionsistrue.
Examples
SIR model:
using JumpProcesses
β = 0.1 / 1000.0;
ν = 0.01;
p = (β, ν)
rate1(u, p, t) = p[1]*u[1]*u[2] # β*S*I
function affect1!(integrator)
integrator.u[1] -= 1 # S -> S - 1
integrator.u[2] += 1 # I -> I + 1
end
jump = ConstantRateJump(rate1, affect1!)
rate2(u, p, t) = p[2]*u[2] # ν*I
function affect2!(integrator)
integrator.u[2] -= 1 # I -> I - 1
integrator.u[3] += 1 # R -> R + 1
end
jump2 = ConstantRateJump(rate2, affect2!)
u₀ = [999, 1, 0]
tspan = (0.0, 250.0)
prob = DiscreteProblem(u₀, tspan, p)
jump_prob = JumpProblem(prob, Direct(), jump, jump2)
sol = solve(jump_prob, SSAStepper())see the tutorial for details.
JumpProcesses.reset_aggregated_jumps! — Functionreset_aggregated_jumps!(integrator, uprev = nothing; update_jump_params=true)Reset the state of jump processes and associated solvers following a change in parameters or such.
Notes
update_jump_params=truewill recalculate the rates stored within any MassActionJump that was built from the parameter vector. If the parameter vector is unchanged, this can safely be set to false to improve performance.
Jump Types
JumpProcesses.ConstantRateJump — Typestruct ConstantRateJump{F1, F2} <: JumpProcesses.AbstractJumpDefines a jump process with a rate (i.e. hazard, intensity, or propensity) that does not explicitly depend on time. More precisely, one where the rate function is constant between the occurrence of jumps. For detailed examples and usage information, see the
Fields
rate: Functionrate(u,p,t)that returns the jump's current rate.affect!: Functionaffect(integrator)that updates the state for one occurrence of the jump.
Examples
Suppose u[1] gives the number of particles and p[1] the probability per time each particle can decay away. A corresponding ConstantRateJump for this jump process is
rate(u,p,t) = p[1]*u[1]
affect!(integrator) = integrator.u[1] -= 1
crj = ConstantRateJump(rate, affect!)Notice, here that rate changes in time, but is constant between the occurrence of jumps (when u[1] will decrease).
JumpProcesses.MassActionJump — Typestruct MassActionJump{T, S, U, V} <: JumpProcesses.AbstractMassActionJumpOptimized representation for ConstantRateJumps that can be represented in mass action form, offering improved performance within jump algorithms compared to ConstantRateJump. For detailed examples and usage information, see the
Constructors
MassActionJump(reactant_stoich, net_stoich; scale_rates = true, param_idxs = nothing)
Here reactant_stoich denotes the reactant stoichiometry for each reaction and net_stoich the net stoichiometry for each reaction.
Fields
scaled_rates: The (scaled) reaction rate constants.reactant_stoch: The reactant stoichiometry vectors.net_stoch: The net stoichiometry vectors.param_mapper: Parameter mapping functor to identify reaction rate constants with parameters inpvectors.
Keyword Arguments
scale_rates = true, whether to rescale the reaction rate constants according to the stoichiometry.nocopy = false, whether theMassActionJumpcan alias thescaled_ratesandreactant_stochfrom the input. Note, ifscale_rates=truethis will potentially modify both of these.param_idxs = nothing, indexes in the parameter vector,JumpProblem.prob.p, that correspond to each reaction's rate.
See the tutorial and main docs for details.
Examples
An SIR model with S + I --> 2I at rate β as the first reaction and I --> R at rate ν as the second reaction can be encoded by
p = (β=1e-4, ν=.01)
u0 = [999, 1, 0] # (S,I,R)
tspan = (0.0, 250.0)
rateidxs = [1, 2] # i.e. [β,ν]
reactant_stoich = [
[1 => 1, 2 => 1], # 1*S and 1*I
[2 => 1] # 1*I
]
net_stoich = [
[1 => -1, 2 => 1], # -1*S and 1*I
[2 => -1, 3 => 1] # -1*I and 1*R
]
maj = MassActionJump(reactant_stoich, net_stoich; param_idxs=rateidxs)
prob = DiscreteProblem(u0, tspan, p)
jprob = JumpProblem(prob, Direct(), maj)Notes
JumpProcesses.VariableRateJump — Typestruct VariableRateJump{R, F, R2, R3, R4, I, T, T2} <: JumpProcesses.AbstractJumpDefines a jump process with a rate (i.e. hazard, intensity, or propensity) that may explicitly depend on time. More precisely, one where the rate function is allowed to change between the occurrence of jumps. For detailed examples and usage information, see the
Note that two types of VariableRateJumps are currently supported, with different performance characteritistics.
A general
VariableRateJumporVariableRateJumpwill refer to one in which onlyrateandaffectfunctions are specified.- These are the most general in what they can represent, but require the use of an
ODEProblemorSDEProblemwhose underlying timestepper handles their evolution in time (via the callback interface). - This is the least performant jump type in simulations.
- These are the most general in what they can represent, but require the use of an
Bounded
VariableRateJumps require passing the keyword argumentsurateandrateinterval, corresponding to functionsurate(u, p, t)andrateinterval(u, p, t), see below. These must calculate a time window over which the rate function is bounded by a constant. Note that it is ok if the rate bound would be violated within the time interval due to a change inuarising from anotherConstantRateJump,MassActionJumpor boundedVariableRateJumpbeing executed, as the chosen aggregator will then handle recalculating the rate bound and interval. However, if the bound could be violated within the time interval due to a change inuarising from continuous dynamics such as a coupled ODE, SDE, or a generalVariableRateJump, bounds should not be given. This ensures the jump is classified as a generalVariableRateJumpand properly handled. One can also optionally provide a lower bound function,lrate(u, p, t), via thelratekeyword argument. This can lead to increased performance. The validity of the lower bound should hold under the same conditions and rate interval asurate.- Bounded
VariableRateJumps can currently be used in theCoevolveaggregator, and can therefore be efficiently simulated in pure-jumpDiscreteProblems using theSSASteppertime-stepper. - These can be substantially more performant than general
VariableRateJumps without the rate bound functions.
- Bounded
Reemphasizing, the additional user provided functions leveraged by bounded VariableRateJumps, urate(u, p, t), rateinterval(u, p, t), and the optional lrate(u, p, t) require that
- For
sin[t, t + rateinterval(u, p, t)], we have thatlrate(u, p, t) <= rate(u, p, s) <= urate(u, p, t). - It is ok if these bounds would be violated during the time window due to another
ConstantRateJump,MassActionJumpor boundedVariableRateJumpoccurring. However, they must remain valid ifuchanges for any other reason (for example, due to continuous dynamics like ODEs, SDEs, or generalVariableRateJumps).
Fields
rate: Functionrate(u,p,t)that returns the jump's current rate given stateu, parameterspand timet.affect!: Functionaffect!(integrator)that updates the state for one occurrence of the jump givenintegrator.lrate: Optional functionlrate(u, p, t)that computes a lower bound on the rate in the intervalttot + rateinterval(u, p, t)at timetgiven stateuand parametersp. This bound must rigorously hold during the time interval as long as anotherConstantRateJump,MassActionJump, or boundedVariableRateJumphas not been sampled. When using aggregators that support boundedVariableRateJumps, currently onlyCoevolve, providing a lower-bound can lead to improved performance.
urate: Optional functionurate(u, p, t)for generalVariableRateJumps, but is required to define a boundedVariableRateJump, which can be used with supporting aggregators, currently onlyCoevolve, and offers improved computational performance. Computes an upper bound for the rate in the intervalttot + rateinterval(u, p, t)at timetgiven stateuand parametersp. This bound must rigorously hold during the time interval as long as anotherConstantRateJump,MassActionJump, or boundedVariableRateJumphas not been sampled.rateinterval: Optional functionrateinterval(u, p, t)for generalVariableRateJumps, but is required to define a boundedVariableRateJump, which can be used with supporting aggregators, currently onlyCoevolve, and offers improved computational performance. Computes the time interval from timetover which theurateandlratebounds will hold,ttot + rateinterval(u, p, t), given stateuand parametersp. This bound must rigorously hold during the time interval as long as anotherConstantRateJump,MassActionJump, or boundedVariableRateJumphas not been sampled.idxsrootfindinterp_pointssave_positionsabstolreltol
Examples
Suppose u[1] gives the number of particles and t*p[1] the probability per time each particle can decay away. A corresponding VariableRateJump for this jump process is
rate(u,p,t) = t*p[1]*u[1]
affect!(integrator) = integrator.u[1] -= 1
vrj = VariableRateJump(rate, affect!)To define a bounded VariableRateJump that can be used with supporting aggregators such as Coevolve, we must define bounds and a rate interval:
rateinterval(u,p,t) = (1 / p[1]) * 2
rate(u,p,t) = t * p[1] * u[1]
lrate(u, p, t) = rate(u, p, t)
urate(u,p,t) = rate(u, p, t + rateinterval(u,p,t))
affect!(integrator) = integrator.u[1] -= 1
vrj = VariableRateJump(rate, affect!; lrate = lrate, urate = urate,
rateinterval = rateinterval)Notes
- When using an aggregator that supports bounded
VariableRateJumps,DiscreteProblemcan be used. Otherwise,ODEProblemorSDEProblemmust be used. - When not using aggregators that support bounded
VariableRateJumps, or when there are generalVariableRateJumps,integrators store an effective state type that wraps the main state vector. SeeExtendedJumpArrayfor details on using this object. In this case allConstantRateJump,VariableRateJumpand callbackaffect!functions receive an integrator withintegrator.uanExtendedJumpArray. - Salis H., Kaznessis Y., Accurate hybrid stochastic simulation of a system of coupled chemical or biochemical reactions, Journal of Chemical Physics, 122 (5), DOI:10.1063/1.1835951 is used for calculating jump times with
VariableRateJumps within ODE/SDE integrators.
JumpProcesses.RegularJump — Typestruct RegularJump{iip, R, C, MD}Representation for encoding rates and multiple simultaneous jumps for use in τ-leaping type methods.
Constructors
RegularJump(rate, c, numjumps; mark_dist = nothing)
Fields
rate: Functionrate!(rate_vals, u, p, t)that returns the current rates, i.e. intensities or propensities, for all possible jumps inrate_vals.
c: Functionc(du, u, p, t, counts, mark)that executes theith jumpcounts[i]times, saving the output indu[i].
numjumps: Number of jumps in the system.mark_dist: A distribution for marks. Not currently used or supported.
Examples
function rate!(out, u, p, t)
out[1] = (0.1 / 1000.0) * u[1] * u[2]
out[2] = 0.01u[2]
nothing
end
const dc = zeros(3,2)
function c(du, u, p, t, counts, mark)
mul!(du, dc, counts)
nothing
end
rj = RegularJump(rate!, c, 2)
## Notes
- `mark_dist` is not currently used or supported in τ-leaping methods.JumpProcesses.JumpSet — Typestruct JumpSet{T1, T2, T3, T4} <: JumpProcesses.AbstractJumpDefines a collection of jumps that should collectively be included in a simulation.
Fields
variable_jumps: Collection ofVariableRateJumpsconstant_jumps: Collection ofConstantRateJumpsregular_jump: Collection ofRegularJumpsmassaction_jump: Collection ofMassActionJumps
Examples
Here we construct two jumps, store them in a JumpSet, and then simulate the resulting process.
using JumpProcesses, OrdinaryDiffEq
rate1(u,p,t) = p[1]
affect1!(integrator) = (integrator.u[1] += 1)
crj = ConstantRateJump(rate1, affect1!)
rate2(u,p,t) = (t/(1+t))*p[2]*u[1]
affect2!(integrator) = (integrator.u[1] -= 1)
vrj = VariableRateJump(rate2, affect2!)
jset = JumpSet(crj, vrj)
f!(du,u,p,t) = (du .= 0)
u0 = [0.0]
p = (20.0, 2.0)
tspan = (0.0, 200.0)
oprob = ODEProblem(f!, u0, tspan, p)
jprob = JumpProblem(oprob, Direct(), jset)
sol = solve(jprob, Tsit5())Aggregator Types
Aggregators are the underlying algorithms used for sampling ConstantRateJumps, MassActionJumps, and VariableRateJumps.
JumpProcesses.Coevolve — TypeAn improvement of the COEVOLVE algorithm for simulating any compound jump process that evolves through time. This method handles variable intensity rates with user-defined bounds and inter-dependent processes. It reduces to NRM when rates are constant. As opposed to COEVOLVE, this method syncs the thinning procedure with the stepper which allows it to handle dependencies on continuous dynamics.
G. A. Zagatti, S. A. Isaacson, C. Rackauckas, V. Ilin, S.-K. Ng and S. Bressan, Extending JumpProcess.jl for fast point process simulation with time-varying intensities, arXiv. doi:10.48550/arXiv.2306.06992.
M. Farajtabar, Y. Wang, M. Gomez-Rodriguez, S. Li, H. Zha, and L. Song, COEVOLVE: a joint point process model for information diffusion and network evolution, Journal of Machine Learning Research 18(1), 1305–1353 (2017). doi: 10.5555/3122009.3122050.
JumpProcesses.Direct — TypeGillespie's Direct method. ConstantRateJump rates and affects are stored in tuples. Fastest for a small (total) number of ConstantRateJumps or MassActionJumps (~10). For larger numbers of possible jumps, use other methods.
Daniel T. Gillespie, A general method for numerically simulating the stochastic time evolution of coupled chemical reactions, Journal of Computational Physics, 22 (4), 403–434 (1976). doi:10.1016/0021-9991(76)90041-3.
JumpProcesses.DirectCR — TypeThe Composition-Rejection Direct method. Performs best relative to other methods for systems with large numbers of jumps with special structure (for example a linear chain of reactions, or jumps corresponding to particles hopping on a grid or graph).
A. Slepoy, A.P. Thompson and S.J. Plimpton, A constant-time kinetic Monte Carlo algorithm for simulation of large biochemical reaction networks, Journal of Chemical Physics, 128 (20), 205101 (2008). doi:10.1063/1.2919546.
S. Mauch and M. Stalzer, Efficient formulations for exact stochastic simulation of chemical systems, ACM Transactions on Computational Biology and Bioinformatics, 8 (1), 27-35 (2010). doi:10.1109/TCBB.2009.47.
JumpProcesses.FRM — TypeGillespie's First Reaction Method. Should not be used for practical applications due to slow performance relative to all other methods.
Daniel T. Gillespie, A general method for numerically simulating the stochastic time evolution of coupled chemical reactions, Journal of Computational Physics, 22 (4), 403–434 (1976). doi:10.1016/0021-9991(76)90041-3.
JumpProcesses.NRM — TypeThe Next Reaction Method. Can significantly outperform Direct for systems with large numbers of jumps and sparse dependency graphs, but is usually slower than one of DirectCR, RSSA, or RSSACR for such systems.
M. A. Gibson and J. Bruck, Efficient exact stochastic simulation of chemical systems with many species and many channels, Journal of Physical Chemistry A, 104 (9), 1876-1889 (2000). doi:10.1021/jp993732q.
JumpProcesses.RDirect — TypeA rejection-based direct method.
JumpProcesses.RSSA — TypeThe Rejection SSA method. One of the best methods for systems with hundreds to many thousands of jumps (along with RSSACR) and sparse dependency graphs.
V. H. Thanh, C. Priami and R. Zunino, Efficient rejection-based simulation of biochemical reactions with stochastic noise and delays, Journal of Chemical Physics, 141 (13), 134116 (2014). doi:10.1063/1.4896985
V. H. Thanh, R. Zunino and C. Priami, On the rejection-based algorithm for simulation and analysis of large-scale reaction networks, Journal of Chemical Physics, 142 (24), 244106 (2015). doi:10.1063/1.4922923.
JumpProcesses.RSSACR — TypeThe Rejection SSA Composition-Rejection method. Often the best performer for systems with tens of thousands of jumps and sparse dependency graphs.
V. H. Thanh, R. Zunino, and C. Priami, Efficient constant-time complexity algorithm for stochastic simulation of large reaction networks, IEEE/ACM Transactions on Computational Biology and Bioinformatics, 14 (3), 657-667 (2017). doi:10.1109/TCBB.2016.2530066.
JumpProcesses.SortingDirect — TypeThe Sorting Direct method. Often the fastest algorithm for smaller to moderate sized systems (tens of jumps), or systems where a few jumps occur much more frequently than others.
J. M. McCollum, G. D. Peterson, C. D. Cox, M. L. Simpson and N. F. Samatova, The sorting direct method for stochastic simulation of biochemical systems with varying reaction execution behavior, Computational Biology and Chemistry, 30 (1), 39049 (2006). doi:10.1016/j.compbiolchem.2005.10.007.
Private API Functions
JumpProcesses.ExtendedJumpArray — Typestruct ExtendedJumpArray{T3<:Number, T1, T<:AbstractArray{T3<:Number, T1}, T2} <: AbstractArray{T3<:Number, 1}Extended state definition used within integrators when there are VariableRateJumps in a system. For detailed examples and usage information, see the
Fields
u: The current state.jump_u: The current rate (i.e. hazard, intensity, or propensity) values for theVariableRateJumps.
Examples
using JumpProcesses, OrdinaryDiffEq
f(du, u, p, t) = du .= 0
rate(u, p, t) = (1+t)*u[1]*u[2]
# suppose we wish to decrease each of the two variables by one
# when a jump occurs
function affect!(integrator)
# Method 1, direct indexing works like normal
integrator.u[1] -= 1
integrator.u[2] -= 1
# Method 2, if we want to broadcast or use array operations we need
# to access integrator.u.u which is the actual state object.
# So equivalently to above we could have said:
# integrator.u.u .-= 1
end
u0 = [10.0, 10.0]
vrj = VariableRateJump(rate, affect!)
oprob = ODEProblem(f, u0, (0.0, 2.0))
jprob = JumpProblem(oprob, Direct(), vrj)
sol = solve(jprob, Tsit5())Notes
If
ueja isa ExtendedJumpArraywithueja.uof sizeNandueja.jump_uof sizenum_variableratejumpsthen# for 1 <= i <= N ueja[i] == ueja.u[i] # for N < i <= (N+num_variableratejumps) ueja[i] == ueja.jump_u[i]In a system with
VariableRateJumps all callback,ConstantRateJump, andVariableRateJumpaffect!functions will receive integrators withintegrator.uanExtendedJumpArray.As such,
affect!functions that wish to modify the state via vector operations should useueja.u.uto obtain the aliased state object.
JumpProcesses.SSAIntegrator — Typemutable struct SSAIntegrator{F, uType, tType, tdirType, P, S, CB, SA, OPT, TS} <: JumpProcesses.AbstractSSAIntegrator{SSAStepper, Nothing, uType, tType}Solution objects for pure jump problems solved via SSAStepper.
Fields
f: The underlyingprob.ffunction. Not currently used.u: The current solution values.t: The current solution time.tprev: The previous time a jump occurred.tdir: The direction time is changing in (must be positive, indicating time is increasing)p: The current parameters.sol: The current solution object.itstop: The next jump time.cb: The jump aggregator callback.saveat: Times to save the solution at.save_everystep: Whether to save every time a jump occurs.save_end: Whether to save at the final step.cur_saveat: Index of the nextsaveattime.opts: Tuple storing callbacks.tstops: User supplied times to step to, useful with callbacks.tstops_idxu_modifiedkeep_steppingalias_tstops: If true, will write tstops into the user-passed arraycopied_tstops: If true indicates we have already allocated the tstops array