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 SBML.Model
:
julia> using SBML
julia> mdl = readSBML("Ec_core_flux1.xml")
SBML.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
⋮ ⋮ ⋱ ⋮ ⋮
Function reference
Helper types
SBML.Maybe
— TypeMaybe{X}
Type shortcut for "X
or nothing" or "nullable X
" in javaspeak. Name got inspired by our functional friends.
SBML.VPtr
— TypeVPtr
A convenience wrapper for "any" (C void
) pointer.
Data structures
SBML.Compartment
— TypeSBML Compartment with sizing information.
SBML.FunctionDefinition
— TypeCustom function definition.
SBML.GPAAnd
— TypeBoolean binary "and" in the association expression
SBML.GPAOr
— TypeBoolean binary "or" in the association expression
SBML.GPARef
— TypeGene product reference in the association expression
SBML.GeneProduct
— TypeGene product metadata.
SBML.GeneProductAssociation
— TypeAbstract type for all kinds of gene product associations
SBML.Math
— TypeA simplified representation of MathML-specified math AST
SBML.MathApply
— TypeFunction application ("call by name", no tricks allowed) in mathematical expression
SBML.MathConst
— TypeA constant identified by name (usually something like pi
, e
or true
) in mathematical expression
SBML.MathIdent
— TypeAn identifier (usually a variable name) in mathematical expression
SBML.MathLambda
— TypeFunction definition (aka "lambda") in mathematical expression
SBML.MathTime
— TypeA special value representing the current time of the simulation, with a special name.
SBML.MathVal
— TypeA literal value (usually a numeric constant) in mathematical expression
SBML.Model
— TypeStructure that collects the model-related data. Contains parameters
, units
, compartments
, species
and reactions
and gene_products
, and additional notes
and annotation
(also present internally in some of the data fields). 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
). Also may contains notes
and annotation
.
SBML.Species
— TypeSpecies metadata – contains a human-readable name
, a compartment
identifier, formula
, charge
, and additional notes
and annotation
.
Base functions
SBML.SBMLVersion
— Methodfunction SBMLVersion()
Get the version of the used SBML library in Julia version format.
SBML.extractModel
— Method" function extractModel(mdl::VPtr)::SBML.Model
Take the SBMLModel_t
pointer and extract all information required to make a valid SBML.Model
structure.
SBML.getAssociation
— Methodfunction getAssociation(x::VPtr)::GeneProductAssociation
Convert a pointer to SBML FbcAssociation_t
to the GeneProductAssociation
tree structure.
SBML.get_optional_bool
— Methodget_optional_bool(x::VPtr, is_sym, get_sym)::Maybe{Bool}
Helper for getting out boolean flags.
SBML.get_optional_double
— Methodget_optional_double(x::VPtr, is_sym, get_sym)::Maybe{Float64}
Helper for getting out C doubles aka Float64s.
SBML.get_optional_int
— Methodget_optional_int(x::VPtr, is_sym, get_sym)::Maybe{UInt}
Helper for getting out unsigned integers.
SBML.get_optional_string
— Methodget_optional_string(x::VPtr, fn_sym)::Maybe{String}
Like get_string
, but returns nothing
instead of throwing an exception.
This is used to get notes and annotations and several other things (see get_notes
, get_annotations
)
SBML.get_string
— Methodget_string(x::VPtr, fn_sym)::Maybe{String}
C-call the SBML function fn_sym
with a single parameter x
, interpret the result as a string and return it, or throw exception in case the pointer is NULL.
SBML.readSBML
— FunctionreadSBML(
fn::String,
sbml_conversion = document -> nothing;
report_severities = ["Fatal", "Error"],
)::SBML.Model
Read the SBML from a XML file in fn
and return the contained SBML.Model
.
The sbml_conversion
is a function that does an in-place modification of the single parameter, which is the C pointer to the loaded SBML document (C type SBMLDocument*
). Several functions for doing that are prepared, including set_level_and_version
, libsbml_convert
, and convert_simplify_math
.
report_severities
switches on and off reporting of certain errors; see the documentation of get_error_messages
for details.
Example
m = readSBML("my_model.xml", doc -> begin
set_level_and_version(3, 1)(doc)
convert_simplify_math(doc)
end)
libsbml
representation converters
The converters are intended to be used as parameters of readSBML
.
SBML.convert_simplify_math
— Functionconvert_simplify_math
Shortcut for libsbml_convert
that expands functions, local parameters, and initial assignments in the SBML document.
SBML.libsbml_convert
— Methodlibsbml_convert(conversion_options::Vector{Pair{String, Dict{String, String}}})
A converter that runs the SBML conversion routine, with specified conversion options. The argument is a vector of pairs to allow specifying the order of conversions.
SBML.libsbml_convert
— Methodlibsbml_convert(converter::String; kwargs...)
Quickly construct a single run of a libsbml
converter from keyword arguments.
Example
readSBML("example.xml", libsbml_convert("stripPackage", package="layout"))
SBML.set_level_and_version
— Methodset_level_and_version(level, version)
A converter to pass into readSBML
that enforces certain SBML level and version.
Data helpers
SBML.check_errors
— Functioncheck_errors(
success::Integer,
doc::VPtr,
error::Exception,
report_severities = ["Fatal", "Error"],
)
If success is a 0-valued Integer
(a logical false
), then call get_error_messages
to show the error messages reported by SBML in the doc
document and throw the error
if they are more than 1. success
is typically the value returned by an SBML C function operating on doc
which returns a boolean flag to signal a successful operation.
SBML.extensive_kinetic_math
— Methodfunction extensive_kinetic_math(
m::SBML.Model,
formula::SBML.Math;
handle_empty_compartment_size = (id::String) -> throw(
DomainError(
"Non-substance-only-unit reference to species `$id' in an unsized compartment `$(m.species[id].compartment)'",
),
),
)
Convert a SBML math formula
to "extensive" kinetic laws, where the references to species that are marked as not having only substance units are converted from amounts to concentrations.
If the data is missing, you can supply a function that adds them. A common way to handle errors is to assume that unsized compartments have volume 1.0 (of whatever units), you can specify that behavior by supplying handle_empty_compartment_size = _ -> 1.0
.
Handling of units in the conversion process is ignored in this version.
SBML.getLBs
— MethodgetLBs(m::SBML.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
— MethodgetOCs(m::SBML.Model)::Vector{Float64}
Extract the vector of objective coefficients of each reaction.
SBML.getS
— Methodfunction getS(m::SBML.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 SBML.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
— MethodgetUBs(m::SBML.Model)::Vector{Tuple{Float64,String}}
Likewise to getLBs
, extract the upper bounds.
SBML.get_error_messages
— Methodget_error_messages(doc::VPtr, error::Exception, report_severities)
Show the error messages reported by SBML in the doc
document and throw the error
if they are more than 1.
report_severities
switches the reporting of certain error types defined by libsbml; you can choose from ["Fatal", "Error", "Warning", "Informational"]
.
SBML.initial_amounts
— Methodinitial_amounts(m::SBML.Model; convert_concentrations = false)
Return initial amounts for each species as a generator of pairs species_name => initial_amount
; the amount is set to nothing
if not available. If convert_concentrations
is true and there is information about initial concentration available together with compartment size, the result is computed from the species' initial concentration.
In the current version, units of the measurements are completely ignored.
Example
# get the initial amounts as dictionary
Dict(initial_amounts(model, convert_concentrations = true))
# remove the empty entries
Dict(k => v for (k,v) in initial_amounts(model) if !isnothing(v))
SBML.initial_concentrations
— Methodinitial_concentrations(m::SBML.Model; convert_amounts = false)
Return initial concentrations of the species in the model. Refer to work-alike initial_amounts
for details.
Math and Symbolics.jl
compatibility
SBML.default_symbolics_constants
— Constantconst default_symbolics_constants::Dict{String, Any}
A dictionary of default constants filled in place of SBML Math constants in the symbolics conversion.
SBML.default_symbolics_mapping
— Constantdefault_symbolics_mapping :: Dict{String,Any}
Default mapping of SBML function names to Julia functions, represented as a dictionary from Strings (SBML names) to anything eval
uable as Julia&Symbolics functions, such as symbols and expressions.
The default mapping only contains the basic SBML functions that are unambiguously represented in Julia; it is supposed to be extended by the user if more functions need to be supported.
Base.convert
— MethodBase.convert(
::Type{Num},
x::SBML.Math;
mapping = default_symbolics_mapping,
convert_time = (x::SBML.MathTime) -> Num(Variable(Symbol(x.id))).val,
convert_const = (x::SBML.MathConst) -> Num(default_symbolics_constants[x.id]),
)
Convert SBML.Math
to Num
type from Symbolics package. The conversion of functions can be customized by supplying a custom mapping; if nothing is supplied, default_symbolics_mapping
that translates basic functions to their Julia equivalents is assumed.
Translation of MathLambda
is not supported by Symbolics.
MathTime
is handled specially, the function from the argument convert_time
is called to possibly specify any desired behavior. By default, it just creates a variable with the same name as the time variable name stored in SBML.
Internal math helpers
SBML.ast_is
— Methodast_is(ast::VPtr, what::Symbol)::Bool
Helper for quickly recognizing kinds of ASTs
SBML.parse_math
— Methodparse_math(ast::VPtr)::Math
This attempts to parse out a decent Julia-esque (Math
AST from a pointer to ASTNode_t
.
SBML.parse_math_children
— Methodparse_math_children(ast::VPtr)::Vector{Math}
Recursively parse all children of an AST node.