Skip to content

Tutorial

This tutorial demonstrates how to import an SBML model with SBMLImporter and simulate it as a JumpProblem (Gillespie/SSA), SDEProblem (chemical Langevin), or ODEProblem (deterministic ODEs). It also describes how to modify parameters and initial conditions, as we as how to index solution trajectories by SBML IDs.

Input: a valid SBML file

SBMLImporter only requires one input: a valid SBML file. SBML files can be created from various sources, such as the graphical interface of COPASI, by converting rule-based BioNetGen models, or by using any other SBML exporting tools. Additionally, a large collection of published SBML models are hosted on BioModels.

This tutorial will use the Brusselator model, whose SBML file can be downloaded from this link.

Importing a model

SBML models are imported with load_SBML:

julia
using SBMLImporter
rn, cb = load_SBML(path_sbml; massaction=true)

The keyword massaction=true indicates that reactions follow chemical mass-action kinetics. This enables efficient jump (SSA-type) simulation.

load_SBML returns two outputs: a Catalyst ReactionSystem (rn) and a CallbackSet (cb). The CallbackSet contains SBML events and callbacks generated from any SBML piecewise expressions. The ReactionSystem contains a map of species initial values (u0) and a map of parameter values (ps), which can be accessed as:

julia
u0 = get_u0_map(rn)
ps = get_parameter_map(rn)

Many modeling tasks can be performed with a Catalyst ReactionSystem. For example, reactions can be inspected with:

julia
using Catalyst
reactions(rn)
CACXY+2XC123XXBCY

INFO

The massaction keyword only affects jump simulations. For SDEProblem and ODEProblem workflows, it can be ignored.

Jump simulations

Jump simulations (e.g. Gillespie/SSA) are run by constructing a JumpProblem:

julia
using JumpProcesses
tspan = (0.0, 10.0)
jprob = JumpProblem(rn, u0, tspan, ps)

The resulting JumpProblem can be simulated with solvers from JumpProcesses.jl, for example SSAStepper:

julia
using Plots
sol = solve(jprob, SSAStepper(), callback=cb)
plot(sol; xlabel="Time [s]", ylabel="Species amount")

Further details are available in the JumpProcesses.jl documentation. A guide to choosing a jump simulation algorithm can be found here.

SDE simulations

Chemical Langevin simulations are run by constructing an SDEProblem from the reaction system:

julia
using StochasticDiffEq
tspan = (0.0, 10.0)
sprob = SDEProblem(rn, u0, tspan, ps)

The SDEProblem can be simulated with solvers from StochasticDiffEq.jl, for example LambaEM:

julia
sol = solve(sprob, LambaEM(), callback=cb)
plot(sol; xlabel="Time [s]", ylabel="Species amount")

For more information, see the StochasticDiffEq.jl documentation.

ODE simulations

Deterministic simulations can be obtained by converting the reaction system into an ODESystem:

julia
using ModelingToolkitBase
sys = mtkcompile(ode_model(rn))

Calling mtkcompile inlines assignment rules into the final equations. As a sanity check, the generated ODEs can be inspected:

julia
equations(sys)
d Y(t)dt=B C X(t)(X(t))2 C Y(t)d X(t)dt=A CC X(t)B C X(t)+(X(t))2 C Y(t)

An ODEProblem can then be constructed:

julia
using OrdinaryDiffEq
tspan = (0.0, 10.0)
oprob = ODEProblem(sys, merge(Dict(u0), Dict(ps)), tspan, jac=true)

With jac=true, the Jacobian is generated symbolically, which often improves performance. The ODEProblem can be simulated with solvers from OrdinaryDiffEq.jl, for example Rodas5P:

julia
sol = solve(oprob, Rodas5P(), callback=cb)
plot(sol; xlabel="Time [s]", ylabel="Species amount")

For more information, see the OrdinaryDiffEq.jl documentation. An ODE solver selection guide is available here.

Importing large models

For large models (>1000 species), symbolic Jacobian construction can be time-consuming. If import time dominates, setting jac=false can help.

Indexing parameters, species, and solutions

Problems imported by load_SBML (ODEProblem, SDEProblem, or JumpProblem) support SymbolicIndexingInterface, enabling model variables to be read and updated by SBML IDs.

Parameter values can be modified via prob.ps. For example, parameter A can be set with:

julia
oprob.ps[:A] = 2.0

Initial conditions can be modified by indexing the problem. For example, the initial value of X can be set with:

julia
oprob[:X] = 2.0

After solving, the same symbolic indexing can be used to retrieve trajectories, e.g. for species X:

julia
sol = solve(oprob, Rodas5P(), callback=cb)
sol[:X]
47-element Vector{Float64}:
  2.0
  3.787570124204133
  5.5166734308533
  8.095067081821867
  9.499430337266517
 10.449730413343405
 10.927931715773303
 11.131387019936806
 11.136329116553842
 11.011995418201597

  2.022372239967697
  2.1028681422396422
  2.06004630958567
  1.9818959085087768
  1.960466594688629
  1.9966554738597053
  2.0190286919500213
  1.9995696819688278
  1.9917855255968444

The same solution indexing applies to other SBML variables (e.g. assignment rules). If inline_assignment_rules=true was passed to load_SBML, assignment rules are inlined into the model equations during import and cannot be accessed this way.

Next steps

For downstream modeling tasks supported by Catalyst ReactionSystems (e.g. bifurcation analysis), see the Catalyst documentation. For parameter estimation of SBML models, see the PEtab.jl documentation

If importing an SBML model fails, please consult the FAQ first. Additional import options and keyword arguments are documented in the API.