API
SBMLImporter.SBML_to_ODESystem
— FunctionSBML_to_ODESystem(path_SBML::AbstractString;
ifelse_to_callback::Bool=true,
write_to_file::Bool=false,
verbose::Bool=true,
return_all::Bool=false,
model_as_string::Bool=false)
Parse an SBML model into a ModelingToolkit ODESystem
and potentially convert events/piecewise to callbacks.
By default, structurally_simplified
is called on the ODESystem
before it is returned.
For information on simulating the ODESystem
, refer to the documentation.
For testing path_SBML
can be the model as a string if model_as_string=true
The number of returned arguments depends on whether the SBML model has events and/or piecewise expressions (see below).
Arguments
path_SBML
: File path to a valid SBML file (level 2 or higher).ifelse_to_callback=true
: Whether to rewriteifelse
(piecewise) expressions to callbacks; recommended for performance.write_to_file=false
: Whether to write the parsed SBML model to a Julia file in the same directory as the SBML file.verbose=true
: Whether or not to display information on the number of return arguments.return_all=true
: Whether or not to return all possible arguments (see below), regardless of whether the model has events.model_as_string=false
: Whether or not the model (path_SBML
) is provided as str, mainly for testing.
Returns
ode_system
: A ModelingToolkitODESystem
that can be converted into anODEProblem
and solved.specie_map
: A species map setting initial values; together with theODESystem
, it can be converted into anODEProblem
.parameter_map
A parameter map setting parameter values; together with theODESystem
, it can be converted into anODEProblem
.cbset
- only for models with events/piecewise expressions: Callbackset (events) for the model.get_tstops
- Only for models with events/piecewise expressions: Function computing time stops for discrete callbacks in thecbset
.ifelse_t0
- Only for models with time-dependent ifelse (piecewise) expressions: Functions checking and adjusting for callback-rewritten piecewise expressions that are active att=t0
.
Examples
# Import and simulate model without events
using SBMLImporter
sys, specie_map, parameter_map = SBML_to_ODESystem(path_SBML)
using OrdinaryDiffEq
tspan = (0, 10.0)
prob = ODEProblem(sys, specie_map, tspan, parameter_map, jac=true)
# Solve ODE with Rodas5P solver
sol = solve(prob, Rodas5P())
# Import a model with events
using SBMLImporter
sys, specie_map, parameter_map, cb, get_tstops = SBML_to_ODESystem(path_SBML)
using OrdinaryDiffEq
tspan = (0, 10.0)
prob = ODEProblem(sys, specie_map, tspan, parameter_map, jac=true)
# Compute event times
tstops = get_tstops(prob.u0, prob.p)
sol = solve(prob, Rodas5P(), tstops=tstops, callback=callbacks)