SBML.jl — load systems biology models from SBML files
This package provides a straightforward way to load model- and simulation-relevant information from SBML files.
The library provides a single function readSBML to load a Model:
julia> using SBML
julia> mdl = readSBML("Ec_core_flux1.xml")
Model(…)
julia> mdl.compartments
2-element Array{String,1}:
"Extra_organism"
"Cytosol"There are several functions to help you with using the data in the usual COBRA-style workflows, such as getS:
julia> metabolites, reactions, S = getS(mdl)
julia> metabolites
77-element Array{String,1}:
"M_succoa_c"
"M_ac_c"
"M_etoh_c"
⋮
julia> S
77×77 SparseArrays.SparseMatrixCSC{Float64,Int64} with 308 stored entries:
[60, 1] = -1.0
[68, 1] = 1.0
[1 , 2] = 1.0
[6 , 2] = -1.0
⋮
[23, 76] = 1.0
[56, 76] = -1.0
[30, 77] = -1.0
[48, 77] = 1.0
julia> Matrix(S)
77×77 Array{Float64,2}:
0.0 1.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 1.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 -1.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 -1.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 -1.0 0.0 0.0 0.0 0.0 0.0
⋮ ⋮ ⋱ ⋮ ⋮ Functions
Data structures
SBML.Model — TypeStructure that collects the model-related data. Contains units, compartments, species and reactions. The contained dictionaries are indexed by identifiers of the corresponding objects.
SBML.Reaction — TypeReaction with stoichiometry that assigns reactants and products their relative consumption/production rates (accessible in field stoichiometry), lower/upper bounds (in tuples lb and ub, with unit names), and objective coefficient (oc).
SBML.Species — TypeSpecies metadata – contains a human-readable name, and a compartment identifier
SBML.UnitPart — TypePart of a measurement unit definition that corresponds to the SBML definition of Unit. For example, the unit "per square megahour", Mh^(-2), is written as:
UnitPart("second", # base SI unit, this says we are measuring time
-2, # exponent, says "per square"
6, # log-10 scale of the unit, says "mega"
1/3600) # second-to-hour multiplierCompound units (such as "volt-amperes" and "dozens of yards per ounce") are built from multiple UnitParts; see the definition of field units in Model.
Base functions
SBML.SBMLVersion — Methodfunction SBMLVersion()Get the version of the used SBML library in Julia version format.
SBML.readSBML — Methodfunction readSBML(fn::String)::ModelRead the SBML from a XML file in fn and return the contained Model.
Data helpers
SBML.getLBs — Methodfunction getLBs(m::Model)::Vector{Tuple{Float64,String}}Extract a vector of lower bounds of reaction rates from the model. All bounds are accompanied with the unit of the corresponding value (this behavior is based on SBML specification).
SBML.getOCs — Methodfunction getOCs(m::Model)::Vector{Float64}Extract the vector of objective coefficients of each reaction.
SBML.getS — Methodfunction getS(m::Model; zeros=spzeros)::Tuple{Vector{String},Vector{String},AbstractMatrix{Float64}}Extract the vector of species (aka metabolite) identifiers, vector of reaction identifiers, and the (dense) stoichiometry matrix from an existing Model. Returns a tuple with these values.
The matrix is sparse by default (initially constructed by SparseArrays.spzeros). You can fill in a custom empty matrix constructed to argument zeros; e.g. running with zeros=zeros will produce a dense matrix.
SBML.getUBs — Methodfunction getUBs(m::Model)::Vector{Tuple{Float64,String}}Likewise to getLBs, extract the upper bounds.