API
SBMLImporter.load_SBML
— Functionload_SBML(path; massaction=false, kwargs...)
Import an SBML model into a ParsedReactionNetwork
that can be simulated as a JumpProblem
for Gillespie simulations, a SDEProblem
for chemical Langevin simulations, or an ODEProblem
for deterministic simulations.
The keyword massaction
should be set to true
if the model reactions follow chemical mass action kinetics. This is because the most efficient JumpProblem
(Gillespie) simulators require mass action jumps, for more details see here. For how to check if a model follows mass action kinetics, see the FAQ in the documentation.
The massaction
keyword argument is only relevant for jump simulations. If you import the model as a SDEProblem
or an ODEProblem
, you can (and should) ignore this keyword argument.
Keyword arguments
complete::Bool=true
: Whether or not to mark the returned CatalystReactionSystem
as complete. Only a complete system can be used for simulations, while only an incomplete system can be composed with other systems. For more details, see the Catalyst.jl documentation.ifelse_to_callback::Bool=true
: Rewriteifelse
(SBML piecewise) expressions to callbacks. This improves simulation runtime and stability. Strongly recomended to set totrue
.inline_assignment_rules::Bool=false
: Inline SBML assignment rules into model equations. Recommended for importing large models. However, if set totrue
, it is not possible to access assignment rule variables viasol[:var]
.inline_kineticlaw_parameters::Bool=true
: Inline (hardcode) SBML kinetic law parameter values into the model equations. Kinetic law parameters are parameters defined only in the kinetic law field of an SBML reaction. It is recommended to inline them as there might be duplicate IDs, e.g., the same kinetic law parameter redefined for different reactions. If they are not inlined, they are promoted to model parameters (with the same status as parameters in the SBML file).write_to_file=false
: Write the parsed SBML model to a Julia file in the same directory as the SBML file.model_as_string::Bool=false
: Whether the model (path
) is provided as a string.
Returns
prn
: AParsedReactionNetwork
that can be converted into aJumpProblem
, aSDEProblem
, or anODEProblem
. Important fields are:prn.rn
: A Catalyst.jlReactionNetwork
.prn.u0
: A value map for setting initial values for model species.prn.p
: A value map for setting model parameter values.
cbset
: ACallbackSet
with SBML events and SBML piecewise functions.
Examples
# Import and simulate model as a JumpProblem
using SBMLImporter, JumpProcesses
prn, cb = load_SBML(path; massaction=true)
tspan = (0.0, 10.0)
dprob = DiscreteProblem(prn.rn, prn.u0, tspan, prn.p)
jprob = JumpProblem(prn.rn, dprob, Direct())
sol = solve(jprob, SSAStepper(), callback=cb)
# Import and simulate model as a SDE
using SBMLImporter, StochasticDiffEq
prn, cb = load_SBML(path)
tspan = (0.0, 10.0)
sprob = SDEProblem(prn.rn, prn.u0, tspan, prn.p)
sol = solve(sprob, LambaEM(), callback=cb)
# Import and simulate model as an ODE
using SBMLImporter, OrdinaryDiffEq
prn, cb = load_SBML(path)
oprob = ODEProblem(prn.rn, prn.u0, tspan, prn.p)
sol = solve(oprob, Rodas5P(), callback=cb)
SBMLImporter.getcompartment
— Functiongetcompartment(s)
For an SBML specie s
, get its associated SBML compartment
.
If s
does not have a compartment, or is not a specie, nothing
is returned.
Example
# Get compartment for first model specie
using SBMLImporter, Catalyst
prn, cb = load_SBML(path_SBML)
sbml_species = species(prn.rn)
c = getcompartment(sbml_species[1])