Data types
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 — TypeVPtrA convenience wrapper for "any" (C void) pointer.
Model data structures
SBML.AlgebraicRule — Typestruct AlgebraicRule <: SBML.RuleSBML algebraic rule.
Fields
math::SBML.Math
SBML.AssignmentRule — Typestruct AssignmentRule <: SBML.RuleSBML assignment rule.
Fields
variable::Stringmath::SBML.Math
SBML.CVTerm — Typestruct CVTerm <: SBML.SBMLObjectRepresentation of a SBML CVTerm, usually carrying Model or Biological qualifier, a list of resources, and possibly nested CV terms.
Fields
biological_qualifier::Union{Nothing, Symbol}model_qualifier::Union{Nothing, Symbol}resource_uris::Vector{String}nested_cvterms::Vector{SBML.CVTerm}
SBML.Compartment — Typestruct Compartment <: SBML.SBMLObjectSBML Compartment with sizing information.
Fields
name::Union{Nothing, String}constant::Union{Nothing, Bool}spatial_dimensions::Union{Nothing, Int64}size::Union{Nothing, Float64}units::Union{Nothing, String}metaid::Union{Nothing, String}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.Constraint — Typestruct Constraint <: SBML.SBMLObjectSBML constraint.
Fields
math::SBML.Mathmessage::String
SBML.Event — Typestruct Event <: SBML.SBMLObjectFields
use_values_from_trigger_time::Boolname::Union{Nothing, String}trigger::Union{Nothing, SBML.Trigger}event_assignments::Union{Nothing, Vector{SBML.EventAssignment}}
SBML.EventAssignment — Typestruct EventAssignment <: SBML.SBMLObjectFields
variable::Stringmath::Union{Nothing, SBML.Math}
SBML.FunctionDefinition — Typestruct FunctionDefinition <: SBML.SBMLObjectCustom function definition.
Fields
name::Union{Nothing, String}metaid::Union{Nothing, String}body::Union{Nothing, SBML.Math}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.GPAAnd — Typestruct GPAAnd <: SBML.GeneProductAssociationBoolean binary "and" in the association expression
Fields
terms::Vector{SBML.GeneProductAssociation}
SBML.GPAOr — Typestruct GPAOr <: SBML.GeneProductAssociationBoolean binary "or" in the association expression
Fields
terms::Vector{SBML.GeneProductAssociation}
SBML.GPARef — Typestruct GPARef <: SBML.GeneProductAssociationGene product reference in the association expression
Fields
gene_product::String
SBML.GeneProduct — Typestruct GeneProduct <: SBML.SBMLObjectGene product metadata.
Fields
label::Stringname::Union{Nothing, String}metaid::Union{Nothing, String}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.GeneProductAssociation — Typeabstract type GeneProductAssociation <: SBML.SBMLObjectAbstract type for all kinds of gene product associations
SBML.Group — Typestruct Group <: SBML.SBMLObjectFields
metaid::Union{Nothing, String}kind::Union{Nothing, String}name::Union{Nothing, String}members::Vector{SBML.Member}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.Math — Typeabstract type Math <: SBML.SBMLObjectA simplified representation of MathML-specified math AST
SBML.MathApply — Typestruct MathApply <: SBML.MathFunction application ("call by name", no tricks allowed) in mathematical expression
Fields
fn::Stringargs::Vector{SBML.Math}
SBML.MathAvogadro — Typestruct MathAvogadro <: SBML.MathA special value representing the Avogadro constant (which is a special named value in SBML).
Fields
id::String
SBML.MathConst — Typestruct MathConst <: SBML.MathA constant identified by name (usually something like pi, e or true) in mathematical expression
Fields
id::String
SBML.MathIdent — Typestruct MathIdent <: SBML.MathAn identifier (usually a variable name) in mathematical expression
Fields
id::String
SBML.MathLambda — Typestruct MathLambda <: SBML.MathFunction definition (aka "lambda") in mathematical expression
Fields
args::Vector{String}body::SBML.Math
SBML.MathTime — Typestruct MathTime <: SBML.MathA special value representing the current time of the simulation, with a special name.
Fields
id::String
SBML.MathVal — Typestruct MathVal{T} <: SBML.MathA literal value (usually a numeric constant) in mathematical expression
Fields
val::Any
SBML.Member — Typestruct Member <: SBML.SBMLObjectFields
id::Union{Nothing, String}metaid::Union{Nothing, String}name::Union{Nothing, String}id_ref::Union{Nothing, String}metaid_ref::Union{Nothing, String}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.Model — Typestruct Model <: SBML.SBMLObjectJulia representation of SBML Model structure, with the reactions, species, units, compartments, and many other things.
Where available, all objects are contained in dictionaries indexed by SBML identifiers.
Fields
parameters::Dict{String, SBML.Parameter}units::Dict{String, SBML.UnitDefinition}compartments::Dict{String, SBML.Compartment}species::Dict{String, SBML.Species}initial_assignments::Dict{String, SBML.Math}rules::Vector{SBML.Rule}constraints::Vector{SBML.Constraint}reactions::Dict{String, SBML.Reaction}objectives::Dict{String, SBML.Objective}active_objective::Union{Nothing, String}gene_products::Dict{String, SBML.GeneProduct}function_definitions::Dict{String, SBML.FunctionDefinition}events::Vector{Pair{Union{Nothing, String}, SBML.Event}}groups::Dict{String, SBML.Group}name::Union{Nothing, String}id::Union{Nothing, String}metaid::Union{Nothing, String}conversion_factor::Union{Nothing, String}area_units::Union{Nothing, String}extent_units::Union{Nothing, String}length_units::Union{Nothing, String}substance_units::Union{Nothing, String}time_units::Union{Nothing, String}volume_units::Union{Nothing, String}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.Objective — Typestruct Objective <: SBML.SBMLObjectFields
type::Stringflux_objectives::Dict{String, Float64}
SBML.Parameter — Typestruct Parameter <: SBML.SBMLObjectRepresentation of SBML Parameter structure, holding a value annotated with units and constantness information.
Fields
name::Union{Nothing, String}value::Union{Nothing, Float64}units::Union{Nothing, String}constant::Union{Nothing, Bool}metaid::Union{Nothing, String}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.RateRule — Typestruct RateRule <: SBML.RuleSBML rate rule.
Fields
variable::Stringmath::SBML.Math
SBML.Reaction — Typestruct Reaction <: SBML.SBMLObjectReaction with stoichiometry that assigns reactants and products their relative consumption/production rates, lower/upper bounds (in tuples lb and ub, with unit names), and objective coefficient (oc). Also may contains notes and annotation.
Fields
name::Union{Nothing, String}reactants::Vector{SBML.SpeciesReference}products::Vector{SBML.SpeciesReference}kinetic_parameters::Dict{String, SBML.Parameter}lower_bound::Union{Nothing, String}upper_bound::Union{Nothing, String}gene_product_association::Union{Nothing, SBML.GeneProductAssociation}kinetic_math::Union{Nothing, SBML.Math}reversible::Boolmetaid::Union{Nothing, String}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.Rule — Typeabstract type Rule <: SBML.SBMLObjectAbstract type representing SBML rules.
SBML.SBMLObject — Typeabstract type SBMLObjectCommon supertype for all SBML.jl objects.
SBML.Species — Typestruct Species <: SBML.SBMLObjectSpecies metadata – contains a human-readable name, a compartment identifier, formula, charge, and additional notes and annotation.
Fields
name::Union{Nothing, String}compartment::Stringboundary_condition::Union{Nothing, Bool}formula::Union{Nothing, String}charge::Union{Nothing, Int64}initial_amount::Union{Nothing, Float64}initial_concentration::Union{Nothing, Float64}substance_units::Union{Nothing, String}conversion_factor::Union{Nothing, String}only_substance_units::Union{Nothing, Bool}constant::Union{Nothing, Bool}metaid::Union{Nothing, String}notes::Union{Nothing, String}annotation::Union{Nothing, String}sbo::Union{Nothing, String}cv_terms::Vector{SBML.CVTerm}
SBML.SpeciesReference — Typestruct SpeciesReference <: SBML.SBMLObjectSBML SpeciesReference.
Fields
id::Union{Nothing, String}species::Stringstoichiometry::Union{Nothing, Float64}constant::Union{Nothing, Bool}
SBML.Trigger — Typestruct Trigger <: SBML.SBMLObjectFields
persistent::Boolinitial_value::Boolmath::Union{Nothing, SBML.Math}
SBML.UnitDefinition — Typestruct UnitDefinition <: SBML.SBMLObjectRepresentation of SBML unit definition, holding the name of the unit and a vector of SBML.UnitParts. See the definition of field units in SBML.Model.
Fields
name::Union{Nothing, String}unit_parts::Vector{SBML.UnitPart}
SBML.UnitPart — Typestruct UnitPart <: SBML.SBMLObjectPart 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:
SBML.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 also SBML.UnitDefinition.
Fields
kind::Stringexponent::Int64scale::Int64multiplier::Float64
Base functions
SBML.SBML — ModuleSBML.jl
| Build status | Documentation |
|---|---|
This is a simple wrap of some of the libSBML functionality, mainly the model loading for purposes of COBRA analysis methods and exploration of ODE system and reaction dynamics.
You might like to try the packages that use SBML.jl; these now include:
- COBREXA.jl, the exascale-ready constraint-based analysis and reconstruction toolkit for finding and modeling steady metabolic fluxes with the models
- Connections to the SciML ModelingToolkit ecosystem:
- SBMLToolkit.jl, for working with the models' reaction dynamics as ODE systems.
- SBMLImporter.jl for importing the SBML models into Catalyst.jl, allowing simulations via Catalyst's
JumpProblem,SDEProblem, orODEProblemand fitting models to data via PEtab.jl.
Other functionality will be added as needed. Feel free to submit a PR that increases the loading "coverage".
Acknowledgements
SBML.jl was developed at the Luxembourg Centre for Systems Biomedicine of the University of Luxembourg (lcsb.uni.lu), and the UCL Research Software Development Group (ucl.ac.uk/arc). The development was supported by European Union's Horizon 2020 Programme under PerMedCoE project (permedcoe.eu) agreement no. 951773, and Chan Zuckerberg Initiative (chanzuckerberg.com) under grant 2020-218578 (5022).
<img src="docs/src/assets/unilu.svg" alt="Uni.lu logo" height="64px"> <img src="docs/src/assets/lcsb.svg" alt="LCSB logo" height="64px"> <img src="docs/src/assets/permedcoe.svg" alt="PerMedCoE logo" height="64px"> <img src="docs/src/assets/ucl.svg" alt="UCL logo" height="64px">
Installation
]add SBML # or
using Pkg; Pkg.add("SBML")Usage
using SBML
m = readSBML("myModel.xml")
# m is now a Model structure with:
m.reactions
m.species
m.compartments
...There are several helper functions, for example you can get a nice list of reactions, metabolites and the stoichiometric matrix as follows:
mets, rxns, S = stoichiometry_matrix(m)SBML.sbml — Methodsbml(sym::Symbol) -> Ptr{Nothing}
A shortcut that loads a function symbol from SBML_jll.
Loading, writing and versioning
SBML._readSBML — Method_readSBML(
symbol::Symbol,
fn::String,
sbml_conversion,
report_severities,
throw_severities
) -> SBML.Model
Internal helper for readSBML.
SBML.get_association — Methodget_association(
x::Ptr{Nothing}
) -> Union{SBML.GPAAnd, SBML.GPAOr, SBML.GPARef}
Convert a pointer to SBML FbcAssociation_t to the GeneProductAssociation tree structure.
SBML.get_cv_terms — Methodget_cv_terms(x::Ptr{Nothing}) -> Vector{SBML.CVTerm}
Shortcut for retrieving SBO term IDs (as strings).
SBML.get_model — Methodget_model(mdl::Ptr{Nothing}) -> SBML.Model
Take the SBMLModel_t pointer and extract all information required to make a valid SBML.Model structure.
SBML.get_optional_bool — Methodget_optional_bool(
x::Ptr{Nothing},
is_sym::Symbol,
get_sym::Symbol
) -> Union{Nothing, Bool}
Helper for getting out boolean flags.
SBML.get_optional_double — Methodget_optional_double(
x::Ptr{Nothing},
is_sym::Symbol,
get_sym::Symbol
) -> Union{Nothing, Float64}
Helper for getting out C doubles aka Float64s.
SBML.get_optional_int — Methodget_optional_int(
x::Ptr{Nothing},
is_sym::Symbol,
get_sym::Symbol
) -> Union{Nothing, Int64}
Helper for getting out unsigned integers.
SBML.get_optional_string — Methodget_optional_string(
x::Ptr{Nothing},
fn_test::Symbol,
fn_sym::Symbol
) -> Union{Nothing, String}
Like get_string, but returns nothing instead of throwing an exception. Also returns values only if fn_test returns true.
SBML.get_optional_string — Methodget_optional_string(
x::Ptr{Nothing},
fn_sym::Symbol
) -> Union{Nothing, 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_parameter — Methodget_parameter(
p::Ptr{Nothing}
) -> Pair{String, SBML.Parameter}
Extract the value of SBML Parameter_t.
SBML.get_sbo_term — Methodget_sbo_term(x::Ptr{Nothing}) -> Union{Nothing, String}
Shortcut for retrieving SBO term IDs (as strings).
SBML.get_string — Methodget_string(x::Ptr{Nothing}, fn_sym::Symbol) -> 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.get_string_from_xmlnode — Methodget_string_from_xmlnode(xmlnode::Ptr{Nothing}) -> String
Helper for converting XML that is not represented by SBML structures to String.
SBML.readSBML — FunctionreadSBML(fn::String; ...) -> SBML.Model
readSBML(
fn::String,
sbml_conversion;
report_severities,
throw_severities
) -> 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, convert_simplify_math and convert_promotelocals_expandfuns.
report_severities and throw_severities switch on and off reporting of certain errors; see the documentation of get_error_messages for details.
To read from a string instead of a file, use readSBMLFromString.
Example
m = readSBML("my_model.xml", doc -> begin
set_level_and_version(3, 1)(doc)
convert_simplify_math(doc)
end)SBML.readSBML — MethodreadSBML(f::Function, args...; kwargs...) -> Any
Shortcut for running a function f on the contents of the SBML file, suitable for the do block syntax. Returns whatever f returns.
Example
readSBML("my_model.xml") do m
@info "model stats" length(m.species) length(m.reactions)
endSBML.readSBMLFromString — FunctionreadSBMLFromString(str::AbstractString; ...) -> SBML.Model
readSBMLFromString(
str::AbstractString,
sbml_conversion;
report_severities,
throw_severities
) -> SBML.Model
Read the SBML from the string str and return the contained SBML.Model.
For the other arguments see the docstring of readSBML, which can be used to read from a file instead of a string.
SBML.writeSBML — MethodwriteSBML(
f::Function,
args...;
kwargs...
) -> Union{Nothing, String}
Shortcut for writing a SBML Model created by function f, suitable for the do block syntax.
Example
writeSBML("my_model.xml") do
compartments = ...
species = ...
Model(; name = "my model", compartments, species)
endSBML.writeSBML — MethodwriteSBML(mdl::SBML.Model, filename::String)
Write the SBML structure in mdl to a file filename.
To write the XML to a string, use writeSBML(mdl::Model).
SBML.writeSBML — MethodwriteSBML(mdl::SBML.Model) -> String
Convert the SBML structure in mdl into XML and return it in a string.
To write directly to a file, use writeSBML(mdl::Model, filename::String).
SBML.Version — MethodVersion() -> VersionNumber
Get the version of the used SBML library in Julia version format.
libsbml representation converters
The converters are intended to be used as parameters of readSBML.
SBML.convert_promotelocals_expandfuns — FunctionShortcut for libsbml_convert that expands functions and local parameters in the SBML document.
SBML.convert_simplify_math — FunctionShortcut for libsbml_convert that expands functions, local parameters, and initial assignments in the SBML document.
SBML.libsbml_convert — Functionlibsbml_convert(
conversion_options::AbstractVector{<:Pair{String, <:AbstractDict{String, String}}}
) -> SBML.var"#20#21"{<:AbstractVector{var"#s56"}, Vector{String}, Vector{String}} where var"#s56"<:(Pair{String, <:AbstractDict{String, String}})
libsbml_convert(
conversion_options::AbstractVector{<:Pair{String, <:AbstractDict{String, String}}},
report_severities
) -> SBML.var"#20#21"{var"#s182", _A, Vector{String}} where {var"#s45"<:(Pair{String, <:AbstractDict{String, String}}), var"#s182"<:AbstractVector{var"#s45"}, _A}
libsbml_convert(
conversion_options::AbstractVector{<:Pair{String, <:AbstractDict{String, String}}},
report_severities,
throw_severities
) -> SBML.var"#20#21"{<:AbstractVector{var"#s37"}} where var"#s37"<:(Pair{String, <:AbstractDict{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. report_severities switches on and off reporting of certain errors; see the documentation of get_error_messages for details.
SBML.libsbml_convert — Methodlibsbml_convert(
converter::String;
report_severities,
throw_severities,
kwargs...
) -> SBML.var"#20#21"{Vector{Pair{String, Dict{String, String}}}, Vector{String}, Vector{String}}
Quickly construct a single run of a libsbml converter from keyword arguments. report_severities switches on and off reporting of certain errors; see the documentation of get_error_messages for details.
Example
readSBML("example.xml", libsbml_convert("stripPackage", package="layout"))SBML.set_level_and_version — Functionset_level_and_version(
level,
version
) -> SBML.var"#18#19"{_A, _B, Vector{String}, Vector{String}} where {_A, _B}
set_level_and_version(
level,
version,
report_severities
) -> SBML.var"#18#19"{_A, _B, _C, Vector{String}} where {_A, _B, _C}
set_level_and_version(
level,
version,
report_severities,
throw_severities
) -> SBML.var"#18#19"
A converter to pass into readSBML that enforces certain SBML level and version. report_severities switches on and off reporting of certain errors; see the documentation of get_error_messages for details.
Helper functions
Data accessors
Base.show — Methodshow(io::IO, _::MIME{Symbol("text/plain")}, m::SBML.Model)
Pretty-printer for a SBML model. Avoids flushing too much stuff to terminal by accident.
SBML.check_errors — Functioncheck_errors(
success::Integer,
doc::Ptr{Nothing},
error::Exception
) -> Union{Nothing, Bool}
check_errors(
success::Integer,
doc::Ptr{Nothing},
error::Exception,
report_severities
) -> Union{Nothing, Bool}
check_errors(
success::Integer,
doc::Ptr{Nothing},
error::Exception,
report_severities,
throw_severities
) -> Union{Nothing, Bool}
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 — Methodextensive_kinetic_math(
m::SBML.Model,
formula::SBML.Math
) -> Any
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. Compartment sizes are referenced by compartment identifiers. A compartment with no obvious definition available in the model (as detected by seemsdefined) is either defaulted as size-less (i.e., size is 1.0) in case it does not have spatial dimensions, or reported as erroneous.
SBML.fbc_flux_objective — Methodfbc_flux_objective(
m::SBML.Model,
oid::String
) -> Vector{Float64}
Get the specified FBC maximization objective from a model, as a vector in the same order as keys(m.reactions).
SBML.flux_bounds — Methodflux_bounds(
m::SBML.Model
) -> Tuple{Vector{Tuple{Float64, String}}, Vector{Tuple{Float64, String}}}
Extract the vectors of lower and upper bounds of reaction rates from the model, in the same order as keys(m.reactions). All bounds are accompanied with the unit of the corresponding value (the behavior is based on SBML specification). Missing bounds are represented by negative/positive infinite values with empty-string unit.
SBML.flux_objective — Methodflux_objective(m::SBML.Model) -> Vector{Float64}
Collect a single maximization objective from FBC, and from kinetic parameters if FBC is not available. Fails if there is more than 1 FBC objective.
Provided for simplicity and compatibility with earlier versions of SBML.jl.
SBML.get_compartment_size — Methodget_compartment_size(
m::SBML.Model,
compartment;
default
) -> Union{Nothing, Float64}
A helper for easily getting out a defaulted compartment size.
SBML.get_error_messages — Methodget_error_messages(
doc::Ptr{Nothing},
error::Exception,
report_severities,
throw_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"].
throw_severities selects errors on which the functions throws the error.
SBML.initial_amounts — Methodinitial_amounts(
m::SBML.Model;
convert_concentrations,
compartment_size
) -> Base.Generator{Dict{String, SBML.Species}, SBML.var"#112#116"{Bool, SBML.var"#114#118"{SBML.Model}}}
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.
The units of measurement are ignored in this computation, but one may reconstruct them from substance_units field of Species structure.
Example
# get the initial amounts as dictionary
Dict(SBML.initial_amounts(model, convert_concentrations = true))
# suppose the compartment size is 10.0 if unspecified
collect(SBML.initial_amounts(
model,
convert_concentrations = true,
compartment_size = comp -> SBML.get_compartment_size(model, comp, 10.0),
))
# remove the empty entries
Dict(k => v for (k,v) in SBML.initial_amounts(model) if !isnothing(v))SBML.initial_concentrations — Methodinitial_concentrations(
m::SBML.Model;
convert_amounts,
compartment_size
) -> Base.Generator{Dict{String, SBML.Species}, SBML.var"#121#125"{Bool, SBML.var"#123#127"{SBML.Model}}}
Return initial concentrations of the species in the model. Refer to work-alike initial_amounts for details.
SBML.isfreein — Methodisfreein(id::String, expr::SBML.Math)Determine if id is used and not bound (aka. free) in expr.
SBML.kinetic_flux_objective — Methodkinetic_flux_objective(
m::SBML.Model
) -> Union{Vector{Float64}, Vector{Union{Nothing, Float64}}, Vector{Nothing}}
Get a kinetic-parameter-specified flux objective from the model, as a vector in the same order as keys(m.reactions).
SBML.mayfirst — Methodmayfirst(args...) -> Any
Helper to get the first non-nothing value from the arguments.
SBML.maylift — Methodmaylift(f, args...) -> Any
Helper to lift a function to work on Maybe, returning nothing whenever there's a nothing in args.
SBML.seemsdefined — Methodisboundbyrules(
id::String,
m::SBML.Model
)Determine if an identifier seems defined or used by any Rules in the model.
SBML.stoichiometry_matrix — Methodstoichiometry_matrix(
m::SBML.Model
) -> Tuple{Vector{String}, Vector{String}, SparseArrays.SparseMatrixCSC{Float64, Int64}}
Extract the vector of species (aka metabolite) identifiers, vector of reaction identifiers, and a sparse stoichiometry matrix (of type SparseMatrixCSC from SparseArrays package) from an existing SBML.Model. Returns a 3-tuple with these values.
SBML.test_suite_url — Methodtest_suite_url(case::Int64; ref, level, version) -> String
Get the URL of a SBML test case file in the sbmlteam/sbml-test-suite GitHub repository. The URL can be downloaded using Downloads.download or any other function. Preferably, users should implement a download cache and verify the consistency of the downloaded files before using them.
Units support
SBML.unitful — Methodunitful(
m::SBML.Model,
val::Tuple{Float64, String},
default_unit::Number
) -> Any
Overload of unitful that uses the default_unit if the unit is not found in the model.
Example
julia> SBML.unitful(mdl, (10.0,"firkin"), 90 * u"lb")
990.0 lbSBML.unitful — Methodunitful(
m::SBML.Model,
val::Tuple{Float64, String},
default_unit::String
) -> Any
Overload of unitful that allows specification of the default_unit by string ID.
SBML.unitful — Methodunitful(m::SBML.Model, val::Tuple{Float64, String}) -> Any
Computes a properly unitful value from a value-unit pair stored in the model m.
SBML.unitful — Methodunitful(u::SBML.UnitDefinition) -> Any
Converts an SBML unit definition (i.e., its vector of UnitParts) to a corresponding Unitful unit.
SBML.unitful — Methodunitful(u::SBML.UnitPart) -> Any
Converts a UnitPart to a corresponding Unitful unit.
The conversion is done according to the formula from SBML L3v2 core manual release 2(section 4.4.2).
SBML.unitful — Methodunitful(u::Vector{SBML.UnitPart}) -> Any
Converts an SBML unit (i.e., a vector of UnitParts) to a corresponding Unitful unit.
Math interpretation
SBML.default_constants — ConstantA dictionary of default constants filled in place of SBML Math constants in the function conversion.
SBML.default_function_mapping — ConstantDefault mapping of SBML function names to Julia functions, represented as a dictionary from Strings (SBML names) to functions.
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.
SBML.interpret_math — Methodinterpret_math(
x::SBML.Math;
map_apply,
map_const,
map_ident,
map_lambda,
map_time,
map_avogadro,
map_value
) -> Any
Recursively interpret SBML.Math type. This can be used to relatively easily traverse and evaluate the SBML math, or translate it into any custom representation, such as Expr or the Num of Symbolics.jl (see the SBML test suite for examples).
By default, the function can convert SBML constants, values and function applications, but identifiers, time values and lambdas are not mapped and throw an error. Similarly SBML function rateOf is undefined, users must to supply their own definition of rateOf that uses the correct derivative.
Internal math helpers
SBML.ast_is — Methodast_is(ast::Ptr{Nothing}, what::Symbol) -> Bool
Helper for quickly recognizing kinds of ASTs
SBML.parse_math — Methodparse_math(
ast::Ptr{Nothing}
) -> Union{SBML.MathApply, SBML.MathAvogadro, SBML.MathConst, SBML.MathIdent, SBML.MathLambda, SBML.MathTime, SBML.MathVal}
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::Ptr{Nothing}) -> Vector{SBML.Math}
Recursively parse all children of an AST node.