Data types

Helper types

SBML.MaybeType
Maybe{X}

Type shortcut for "X or nothing" or "nullable X" in javaspeak. Name got inspired by our functional friends.

source
SBML.VPtrType
VPtr

A convenience wrapper for "any" (C void) pointer.

source

Model data structures

SBML.AssignmentRuleType
struct AssignmentRule <: SBML.Rule

SBML assignment rule.

Fields

  • variable::String

  • math::SBML.Math

source
SBML.CVTermType
struct CVTerm

Representation 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}

source
SBML.CompartmentType
struct Compartment

SBML 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}

source
SBML.EventType
struct Event

Fields

  • use_values_from_trigger_time::Int32

  • name::Union{Nothing, String}

  • trigger::Union{Nothing, SBML.Trigger}

  • event_assignments::Union{Nothing, Vector{SBML.EventAssignment}}

source
SBML.FunctionDefinitionType
struct FunctionDefinition

Custom function definition.

Fields

  • name::Union{Nothing, String}

  • body::Union{Nothing, SBML.Math}

  • notes::Union{Nothing, String}

  • annotation::Union{Nothing, String}

source
SBML.GPAAndType
struct GPAAnd <: SBML.GeneProductAssociation

Boolean binary "and" in the association expression

Fields

  • terms::Vector{SBML.GeneProductAssociation}
source
SBML.GPAOrType
struct GPAOr <: SBML.GeneProductAssociation

Boolean binary "or" in the association expression

Fields

  • terms::Vector{SBML.GeneProductAssociation}
source
SBML.GPARefType
struct GPARef <: SBML.GeneProductAssociation

Gene product reference in the association expression

Fields

  • gene_product::String
source
SBML.GeneProductType
struct GeneProduct

Gene product metadata.

Fields

  • label::String

  • name::Union{Nothing, String}

  • metaid::Union{Nothing, String}

  • notes::Union{Nothing, String}

  • annotation::Union{Nothing, String}

  • sbo::Union{Nothing, String}

  • cv_terms::Vector{SBML.CVTerm}

source
SBML.MathType

A simplified representation of MathML-specified math AST

source
SBML.MathApplyType
struct MathApply <: SBML.Math

Function application ("call by name", no tricks allowed) in mathematical expression

Fields

  • fn::String

  • args::Vector{SBML.Math}

source
SBML.MathAvogadroType
struct MathAvogadro <: SBML.Math

A special value representing the Avogadro constant (which is a special named value in SBML).

Fields

  • id::String
source
SBML.MathConstType
struct MathConst <: SBML.Math

A constant identified by name (usually something like pi, e or true) in mathematical expression

Fields

  • id::String
source
SBML.MathIdentType
struct MathIdent <: SBML.Math

An identifier (usually a variable name) in mathematical expression

Fields

  • id::String
source
SBML.MathLambdaType
struct MathLambda <: SBML.Math

Function definition (aka "lambda") in mathematical expression

Fields

  • args::Vector{String}

  • body::SBML.Math

source
SBML.MathTimeType
struct MathTime <: SBML.Math

A special value representing the current time of the simulation, with a special name.

Fields

  • id::String
source
SBML.MathValType
struct MathVal{T} <: SBML.Math

A literal value (usually a numeric constant) in mathematical expression

Fields

  • val::Any
source
SBML.ModelType
struct Model

Structure 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.

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::Dict{String, SBML.Event}

  • 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}

source
SBML.ObjectiveType
struct Objective

Fields

  • type::String

  • flux_objectives::Dict{String, Float64}

source
SBML.ParameterType
struct Parameter

Representation 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}

  • sbo::Union{Nothing, String}

  • cv_terms::Vector{SBML.CVTerm}

source
SBML.RateRuleType
struct RateRule <: SBML.Rule

SBML rate rule.

Fields

  • variable::String

  • math::SBML.Math

source
SBML.ReactionType
struct Reaction

Reaction 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::Bool

  • metaid::Union{Nothing, String}

  • notes::Union{Nothing, String}

  • annotation::Union{Nothing, String}

  • sbo::Union{Nothing, String}

  • cv_terms::Vector{SBML.CVTerm}

source
SBML.RuleType
abstract type Rule

Abstract type representing SBML rules.

source
SBML.SpeciesType
struct Species

Species metadata – contains a human-readable name, a compartment identifier, formula, charge, and additional notes and annotation.

Fields

  • name::Union{Nothing, String}

  • compartment::String

  • boundary_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}

  • 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}

source
SBML.SpeciesReferenceType
struct SpeciesReference

SBML SpeciesReference.

Fields

  • id::Union{Nothing, String}

  • species::String

  • stoichiometry::Union{Nothing, Float64}

  • constant::Union{Nothing, Bool}

source
SBML.TriggerType
struct Trigger

Fields

  • persistent::Bool

  • initial_value::Bool

  • math::Union{Nothing, SBML.Math}

source
SBML.UnitDefinitionType
struct UnitDefinition

Representation 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}

source
SBML.UnitPartType
struct UnitPart

Part 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 multiplier

Compound units (such as "volt-amperes" and "dozens of yards per ounce") are built from multiple UnitParts. See also SBML.UnitDefinition.

Fields

  • kind::String

  • exponent::Int64

  • scale::Int64

  • multiplier::Float64

source

Base functions

SBML.SBMLModule

SBML.jl

Build statusDocumentationStats
CI statusstable documentation dev documentationSBML Downloads

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
  • SBMLToolkit.jl, for working with the reaction dynamics of the models as ODE systems, well connected to the SciML ModelingToolkit ecosystem.

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 (uni.lu/lcsb), 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)
source
SBML.sbmlMethod
sbml(sym::Symbol) -> Ptr{Nothing}

A shortcut that loads a function symbol from SBML_jll.

source

Loading, writing and versioning

SBML._readSBMLMethod
_readSBML(
    symbol::Symbol,
    fn::String,
    sbml_conversion,
    report_severities
) -> SBML.Model

Internal helper for readSBML.

source
SBML.get_associationMethod
get_association(
    x::Ptr{Nothing}
) -> Union{SBML.GPAAnd, SBML.GPAOr, SBML.GPARef}

Convert a pointer to SBML FbcAssociation_t to the GeneProductAssociation tree structure.

source
SBML.get_cv_termsMethod
get_cv_terms(x::Ptr{Nothing}) -> Vector{SBML.CVTerm}

Shortcut for retrieving SBO term IDs (as strings).

source
SBML.get_modelMethod
get_model(mdl::Ptr{Nothing}) -> SBML.Model

Take the SBMLModel_t pointer and extract all information required to make a valid SBML.Model structure.

source
SBML.get_optional_boolMethod
get_optional_bool(
    x::Ptr{Nothing},
    is_sym::Symbol,
    get_sym::Symbol
) -> Union{Nothing, Bool}

Helper for getting out boolean flags.

source
SBML.get_optional_doubleMethod
get_optional_double(
    x::Ptr{Nothing},
    is_sym::Symbol,
    get_sym::Symbol
) -> Union{Nothing, Float64}

Helper for getting out C doubles aka Float64s.

source
SBML.get_optional_intMethod
get_optional_int(
    x::Ptr{Nothing},
    is_sym::Symbol,
    get_sym::Symbol
) -> Union{Nothing, Int64}

Helper for getting out unsigned integers.

source
SBML.get_optional_stringMethod
get_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.

source
SBML.get_optional_stringMethod
get_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)

source
SBML.get_parameterMethod
get_parameter(
    p::Ptr{Nothing}
) -> Pair{String, SBML.Parameter}

Extract the value of SBML Parameter_t.

source
SBML.get_sbo_termMethod
get_sbo_term(x::Ptr{Nothing}) -> Union{Nothing, String}

Shortcut for retrieving SBO term IDs (as strings).

source
SBML.get_stringMethod
get_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.

source
SBML.get_string_from_xmlnodeMethod
get_string_from_xmlnode(xmlnode::Ptr{Nothing}) -> String

Helper for converting XML that is not represented by SBML structures to String.

source
SBML.readSBMLFunction
readSBML(fn::String) -> SBML.Model
readSBML(
    fn::String,
    sbml_conversion;
    report_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, and convert_simplify_math.

report_severities switches 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)
source
SBML.readSBMLFromStringFunction
readSBMLFromString(str::AbstractString) -> SBML.Model
readSBMLFromString(
    str::AbstractString,
    sbml_conversion;
    report_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.

source
SBML.writeSBMLMethod
writeSBML(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).

source
SBML.writeSBMLMethod
writeSBML(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).

source
SBML.VersionMethod
Version() -> VersionNumber

Get the version of the used SBML library in Julia version format.

source

libsbml representation converters

The converters are intended to be used as parameters of readSBML.

SBML.libsbml_convertFunction
libsbml_convert(
    conversion_options::AbstractVector{<:Pair{String, <:AbstractDict{String, String}}}
) -> SBML.var"#18#19"{_A, Vector{String}} where _A
libsbml_convert(
    conversion_options::AbstractVector{<:Pair{String, <:AbstractDict{String, String}}},
    report_severities
) -> SBML.var"#18#19"

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.

source
SBML.libsbml_convertFunction
libsbml_convert(
    converter::String
) -> SBML.var"#18#19"{Vector{Pair{String, Dict{String, String}}}, Vector{String}}
libsbml_convert(
    converter::String,
    report_severities;
    kwargs...
) -> SBML.var"#18#19"{Vector{Pair{String, Dict{String, 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"))
source
SBML.set_level_and_versionFunction
set_level_and_version(
    level,
    version
) -> SBML.var"#16#17"{_A, _B, Vector{String}} where {_A, _B}
set_level_and_version(
    level,
    version,
    report_severities
) -> SBML.var"#16#17"

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.

source

Helper functions

Data accessors

Base.showMethod
show(io::IO, _::MIME{Symbol("text/plain")}, m::SBML.Model)

Pretty-printer for a SBML model. Avoids flushing too much stuff to terminal by accident.

source
SBML.check_errorsFunction
check_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}

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.

source
SBML.extensive_kinetic_mathMethod
extensive_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.

source
SBML.fbc_flux_objectiveMethod
fbc_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).

source
SBML.flux_boundsMethod
flux_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.

source
SBML.flux_objectiveMethod
flux_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.

source
SBML.get_compartment_sizeMethod
get_compartment_size(
    m::SBML.Model,
    compartment;
    default
) -> Union{Nothing, Float64}

A helper for easily getting out a defaulted compartment size.

source
SBML.get_error_messagesMethod
get_error_messages(
    doc::Ptr{Nothing},
    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"].

source
SBML.initial_amountsMethod
initial_amounts(
    m::SBML.Model;
    convert_concentrations,
    compartment_size
) -> Base.Generator{Dict{String, SBML.Species}, SBML.var"#104#108"{Bool, SBML.var"#106#110"{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))
source
SBML.initial_concentrationsMethod
initial_concentrations(
    m::SBML.Model;
    convert_amounts,
    compartment_size
) -> Base.Generator{Dict{String, SBML.Species}, SBML.var"#113#117"{Bool, SBML.var"#115#119"{SBML.Model}}}

Return initial concentrations of the species in the model. Refer to work-alike initial_amounts for details.

source
SBML.isfreeinMethod
isfreein(id::String, expr::SBML.Math)

Determine if id is used and not bound (aka. free) in expr.

source
SBML.kinetic_flux_objectiveMethod
kinetic_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).

source
SBML.mayfirstMethod
mayfirst(args...) -> Any

Helper to get the first non-nothing value from the arguments.

source
SBML.mayliftMethod
maylift(f, args::Union{Nothing, X} where X...) -> Any

Helper to lift a function to work on Maybe, returning nothing whenever there's a nothing in args.

source
SBML.seemsdefinedMethod
isboundbyrules(
    id::String,
    m::SBML.Model
)

Determine if an identifier seems defined or used by any Rules in the model.

source
SBML.stoichiometry_matrixMethod
stoichiometry_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.

source

Units support

SBML.unitfulMethod
unitful(
    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 lb
source
SBML.unitfulMethod
unitful(
    m::SBML.Model,
    val::Tuple{Float64, String},
    default_unit::String
) -> Any

Overload of unitful that allows specification of the default_unit by string ID.

source
SBML.unitfulMethod
unitful(m::SBML.Model, val::Tuple{Float64, String}) -> Any

Computes a properly unitful value from a value-unit pair stored in the model m.

source
SBML.unitfulMethod
unitful(u::SBML.UnitDefinition) -> Any

Converts an SBML unit definition (i.e., its vector of UnitParts) to a corresponding Unitful unit.

source
SBML.unitfulMethod
unitful(u::Vector{SBML.UnitPart}) -> Any

Converts an SBML unit (i.e., a vector of UnitParts) to a corresponding Unitful unit.

source

Math interpretation

SBML.default_constantsConstant

A dictionary of default constants filled in place of SBML Math constants in the function conversion.

source
SBML.default_function_mappingConstant

Default 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.

source
SBML.interpret_mathMethod
interpret_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.

source

Internal math helpers

SBML.ast_isMethod
ast_is(ast::Ptr{Nothing}, what::Symbol) -> Bool

Helper for quickly recognizing kinds of ASTs

source
SBML.parse_mathMethod
parse_math(ast::Ptr{Nothing}) -> Any

This attempts to parse out a decent Julia-esque Math AST from a pointer to ASTNode_t.

source