API reference
SbmlInterface.getinitialconditions
— Methodgetinitialconditions(model)
Create array of initial conditions from a libsbml Model.
Example
julia> u0 = getinitialconditions(model)
2-element Array{Pair{Num,Float64},1}:
A(t) => 1.0
B(t) => 0.0
SbmlInterface.getmodel
— Methodgetmodel(sbmlfile::String)
Given an sbmlfile
return a PyObject (libsbml.Model).
Example
julia> model = getmodel(SBML_FILE)
PyObject <Model conversion_reaction_0 "Conversion Reaction 0">
SbmlInterface.getparameters
— Methodgetparameters(model)
Create array of parameters from a libsbml Model.
Example
julia> pars = getparameters(model)
3-element Array{Pair{Num,Float64},1}:
k1 => 0.8
k2 => 0.6
compartment => 1.0
SbmlInterface.getreactions
— Methodgetreactions(model)::Array
Create array of Reactions from a libsbml model.
Example
julia> rxs = getreactions(model)
2-element Array{ModelingToolkit.Reaction,1}:
ModelingToolkit.Reaction{Any,Float64}(compartment*k1*A(t), SymbolicUtils.Term{Real}[A(t)],
SymbolicUtils.Term{Real}[B(t)], [1.0], [1.0],
Pair{Any,Float64}[B(t) => 1.0, A(t) => -1.0], true)
ModelingToolkit.Reaction{Any,Float64}(compartment*k2*B(t), SymbolicUtils.Term{Real}[B(t)],
SymbolicUtils.Term{Real}[A(t)], [1.0], [1.0],
Pair{Any,Float64}[B(t) => -1.0, A(t) => 1.0], true)
SbmlInterface.sbml2odeproblem
— Methodsbml2odeproblem(sbmlfile::String,tspan;jac::Bool=true)
Convert model from SBML file to an ODEProblem
and return it.
Arguments
sbmlfile::String
: path to SBML file.tspan
: The solution u(t) will be computed for tspan[1] ≤ t ≤ tspan[2].jac::Bool=true
: Jacobian matrix.
Example
julia> prob = sbml2odeproblem(SBML_FILE,(0.0,10.0))
ODEProblem with uType Array{Float64,1} and tType Float64. In-place: true
timespan: (0.0, 10.0)
u0: [0.0, 1.0]
SbmlInterface.sbml2odesystem
— Methodsbml2odesystem(sbmlfile::String)
Given an sbmlfile
return an ODESystem
with default_u0
and default_p
set.
Example
julia> sys = sbml2odesystem(SBML_FILE)
Model ##ODESystem#270 with 2 equations
States (2):
B(t) [defaults to 0.0]
A(t) [defaults to 1.0]
Parameters (3):
k2 [defaults to 0.6]
compartment [defaults to 1.0]
k1 [defaults to 0.8]
SbmlInterface.simulatesbml
— Methodsimulatesbml(sbmlfile::String,tspan;saveat::Real=1.,jac::Bool=true,
solver=OrdinaryDiffEq.Tsit5())
Simulate model from SBML file and return an ODESolution.
Arguments
sbmlfile::String
: path to SBML file.tspan
: The solution u(t) will be computed for tspan[1] ≤ t ≤ tspan[2].saveat::Real=1.
: time intervals of u(t) output.jac::Bool=true
: Jacobian matrix.solver=OrdinaryDiffEq.Tsit5()
: Integration algorithm.
Example
julia> sol = simulatesbml("mymodel.sbml",(0.0,1.0))
retcode: Success
Interpolation: 1st order linear
t: 2-element Array{Float64,1}:
0.0
1.0
u: 2-element Array{Array{Float64,1},1}:
[0.0, 1.0]
[0.4, 0.6]