Types

Base types

COBREXA._defaultMethod
_default(d::T, x::Maybe{T})::T where {T}

Fold the Maybe{T} down to T by defaulting.

source
COBREXA.AnnotationsType
Annotations = Dict{String,Vector{String}}

Dictionary used to store (possible multiple) standardized annotations of something, such as a Metabolite and a Reaction.

Example

Annotations("PubChem" => ["CID12345", "CID54321"])
source
COBREXA.GeneAssociationType
GeneAssociation = Vector{Vector{String}}

An association to genes, represented as a logical formula in a positive disjunctive normal form (DNF). (The 2nd-level vectors of strings are connected by "and" to form conjunctions, and the 1st-level vectors of these conjunctions are connected by "or" to form the DNF.)

source
COBREXA.MetabolicModelType
abstract type MetabolicModel end

A helper supertype that wraps everything usable as a linear-like model for COBREXA functions.

If you want your model type to work with COBREXA, add the MetabolicModel as its supertype, and implement the accessor functions. Accessors reactions, metabolites, stoichiometry, bounds and objective must be implemented; others are not mandatory and default to safe "empty" values.

source
COBREXA.NotesType
Notes = Dict{String,Vector{String}}

Free-form notes about something (e.g. a Gene), categorized by "topic".

source
COBREXA.balanceMethod
balance(a::MetabolicModel)::SparseVec

Get the sparse balance vector of a model (ie. the b from S x = b).

source
COBREXA.boundsMethod
bounds(a::MetabolicModel)::Tuple{Vector{Float64},Vector{Float64}}

Get the lower and upper flux bounds of a model.

source
COBREXA.couplingMethod
coupling(a::MetabolicModel)::SparseMat

Get a matrix of coupling constraint definitions of a model. By default, there is no coupling in the models.

source
COBREXA.coupling_boundsMethod
coupling_bounds(a::MetabolicModel)::Tuple{Vector{Float64},Vector{Float64}}

Get the lower and upper bounds for each coupling bound in a model, as specified by coupling. By default, the model does not have any coupling bounds.

source
COBREXA.gene_annotationsMethod
gene_annotations(a::MetabolicModel, gene_id::String)::Annotations

Return standardized names that identify the corresponding gene or product. The dictionary assigns vectors of possible identifiers to identifier system names, e.g. "PDB" => ["PROT01"].

source
COBREXA.gene_nameMethod
gene_name(model::MetabolicModel, gid::String)

Return the name of gene with ID gid.

source
COBREXA.gene_notesMethod
gene_notes(model::MetabolicModel, gene_id::String)::Notes

Return the notes associated with the gene gene_id in model.

source
COBREXA.genesMethod
genes(a::MetabolicModel)::Vector{String}

Return identifiers of all genes contained in the model. By default, there are no genes.

In SBML, these are usually called "gene products" but we write genes for simplicity.

source
COBREXA.metabolite_annotationsMethod
metabolite_annotations(a::MetabolicModel, metabolite_id::String)::Annotations

Return standardized names that may help to reliably identify the metabolite. The dictionary assigns vectors of possible identifiers to identifier system names, e.g. "ChEMBL" => ["123"] or "PubChem" => ["CID123", "CID654645645"].

source
COBREXA.metabolite_chargeMethod
metabolite_charge(model::MetabolicModel, metabolite_id::String)::Maybe{Int}

Return the charge associated with metabolite metabolite_id in model. Returns nothing if charge not present.

source
COBREXA.metabolite_compartmentMethod
metabolite_compartment(model::MetabolicModel, metabolite_id::String)::Maybe{String}

Return the compartment of metabolite metabolite_id in model if it is assigned. If not, return nothing.

source
COBREXA.metabolite_formulaMethod
metabolite_formula(
    a::MetabolicModel,
    metabolite_id::String,
)::Maybe{MetaboliteFormula}

Return the formula of metabolite metabolite_id in model. Return nothing in case the formula is not known or irrelevant.

source
COBREXA.metabolite_notesMethod
metabolite_notes(model::MetabolicModel, metabolite_id::String)::Notes

Return the notes associated with metabolite reaction_id in model.

source
COBREXA.metabolitesMethod
metabolites(a::MetabolicModel)::Vector{String}

Return a vector of metabolite identifiers in a model.

source
COBREXA.n_genesMethod
n_genes(a::MetabolicModel)::Int

Return the number of genes in the model (as returned by genes). If you just need the number of the genes, this may be much more efficient than calling genes and measuring the array.

source
COBREXA.precache!Method
precache!(a::MetabolicModel)::Nothing

Do whatever is feasible to get the model into a state that can be read from as-quickly-as-possible. This may include e.g. generating helper index structures and loading delayed parts of the model from disk. The model should be modified "transparently" in-place. Analysis functions call this right before applying modifications or converting the model to the optimization model using make_optimization_model; usually on the same machine where the optimizers (and, generally, the core analysis algorithms) will run. The calls are done in a good hope that the performance will be improved.

By default, it should be safe to do nothing.

source
COBREXA.reaction_annotationsMethod
reaction_annotations(a::MetabolicModel, reaction_id::String)::Annotations

Return standardized names that may help identifying the reaction. The dictionary assigns vectors of possible identifiers to identifier system names, e.g. "Reactome" => ["reactomeID123"].

source
COBREXA.reaction_gene_associationMethod
reaction_gene_association(a::MetabolicModel, gene_id::String)::Maybe{GeneAssociation}

Returns the sets of genes that need to be present so that the reaction can work (technically, a DNF on gene availability, with positive atoms only).

For simplicity, nothing may be returned, meaning that the reaction always takes place. (in DNF, that would be equivalent to returning [[]].)

source
COBREXA.reaction_notesMethod
reaction_notes(model::MetabolicModel, reaction_id::String)::Notes

Return the notes associated with reaction reaction_id in model.

source
COBREXA.reaction_stoichiometryMethod
reaction_stoichiometry(model::MetaboliteModel, rid::String)::Dict{String, Float64}

Return the stoichiometry of reaction with ID rid in the model. The dictionary maps the metabolite IDs to their stoichiometric coefficients.

source
COBREXA.reaction_subsystemMethod
reaction_subsystem(model::MetabolicModel, reaction_id::String)::Maybe{String}

Return the subsystem of reaction reaction_id in model if it is assigned. If not, return nothing.

source
COBREXA.reactionsMethod
reactions(a::MetabolicModel)::Vector{String}

Return a vector of reaction identifiers in a model.

source

Model types and contents

COBREXA.CoreModelType
struct CoreModel <: MetabolicModel

A "bare bones" core linear optimization problem of the form, with reaction and metabolite names.

min c^T x
s.t. S x = b
      xₗ ≤ x ≤ xᵤ
source
Base.convertMethod
Base.convert(::Type{CoreModel}, m::M) where {M <: MetabolicModel}

Make a CoreModel out of any compatible model type.

source
COBREXA.boundsMethod
bounds(a::CoreModel)::Tuple{Vector{Float64},Vector{Float64}}

CoreModel flux bounds.

source
COBREXA.CoreModelCoupledType
struct CoreModelCoupled <: MetabolicModel

The linear model with additional coupling constraints in the form

    cₗ ≤ C x ≤ cᵤ
source
Base.convertMethod
Base.convert(::Type{CoreModelCoupled}, mm::MetabolicModel)

Make a CoreModelCoupled out of any compatible model type.

source
COBREXA.couplingMethod
coupling(a::CoreModelCoupled)::SparseMat

Coupling constraint matrix for a CoreModelCoupled.

source
COBREXA.coupling_boundsMethod
coupling_bounds(a::CoreModelCoupled)::Tuple{Vector{Float64},Vector{Float64}}

Coupling bounds for a CoreModelCoupled.

source
COBREXA.FluxSummaryType
FluxSummary

A struct used to store summary information about the solution of a constraint based analysis result.

source
COBREXA.flux_summaryMethod
flux_summary(flux_result::Dict{String, Float64};
    exclude_exchanges = false,
    exchange_prefixes = _constants.exchange_prefixes,
    biomass_strings = _constants.biomass_strings,
    exclude_biomass = false,
    small_flux_bound = 1.0/_constants.default_reaction_bound^2,
    large_flux_bound = _constants.default_reaction_bound,
    keep_unbounded = false,
)::FluxSummary

Summarize a dictionary of fluxes into small, useful representation of the most important information contained. Useful for pretty-printing and quickly exploring the results. Internally this function uses looks_like_biomass_reaction and looks_like_exchange_reaction. The corresponding keyword arguments passed to these functions. Use this if your model has non-standard ids for reactions. Fluxes smaller than small_flux_bound are not stored, while fluxes larger than large_flux_bound are only stored if keep_unbounded is true.

Example

julia> sol = flux_dict(flux_balance_analysis(model, Tulip.Optimizer))
julia> fr = flux_summary(sol)
Biomass:
  BIOMASS_Ecoli_core_w_GAM: 0.8739
Import:
  EX_o2_e:     -21.7995
  EX_glc__D_e: -10.0
  EX_nh4_e:    -4.7653
  EX_pi_e:     -3.2149
Export:
  EX_h_e:      17.5309
  EX_co2_e:    22.8098
  EX_h2o_e:    29.1758
source
COBREXA.flux_variability_summaryMethod
flux_variability_summary(flux_result::Tuple{Dict{String, Dict{String, Float64}}, Dict{String, Dict{String, Float64}}};
    exclude_exchanges = false,
    exchange_prefixes = _constants.exchange_prefixes,
    biomass_strings = _constants.biomass_strings,
    exclude_biomass = false,
    )::FluxVariabilitySummary

Summarize a dictionary of flux dictionaries obtained eg. from flux_variability_analysis_dict. The simplified summary representation is useful for pretty-printing and easily showing the most important results.

Internally this function uses looks_like_biomass_reaction and looks_like_exchange_reaction. The corresponding keyword arguments are passed to these functions. Use this if your model has an uncommon naming of reactions.

Example

julia> sol = flux_variability_analysis_dict(model, Gurobi.Optimizer; bounds = objective_bounds(0.99))
julia> flux_res = flux_variability_summary(sol)
Biomass                     Lower bound   Upper bound
  BIOMASS_Ecoli_core_w_GAM: 0.8652        0.8652
Exchange
  EX_h2o_e:                 28.34         28.34
  EX_co2_e:                 22.0377       22.0377
  EX_o2_e:                  -22.1815      -22.1815
  EX_h_e:                   17.3556       17.3556
  EX_glc__D_e:              -10.0         -10.0
  EX_nh4_e:                 -4.8448       -4.8448
  EX_pi_e:                  -3.2149       -3.2149
  EX_for_e:                 0.0           0.0
  ...                       ...           ...
source
COBREXA.GeneType

Gene struct.

Fields

id :: String
name :: Maybe{String}
notes :: Dict{String, Vector{String}}
annotation :: Dict{String, Union{Vector{String}, String}}
source
COBREXA.JSONModelType
struct JSONModel <: MetabolicModel
    json::Dict{String,Any}
    rxn_index::Dict{String,Int}
    rxns::Vector{Any}
    met_index::Dict{String,Int}
    mets::Vector{Any}
    gene_index::Dict{String,Int}
    genes::Vector{Any}
end

A struct used to store the contents of a JSON model, i.e. a model read from a file ending with .json. These model files typically store all the model data in arrays of JSON objects (represented in Julia as vectors of dictionaries).

Usually, not all of the fields of the input JSON can be easily represented when converting to other models, care should be taken to avoid losing information.

Direct work with the json structure is not very efficient; the model structure therefore caches some of the internal structure in the extra fields. The single-parameter JSONModel constructor creates these caches correctly from the json. The model structure is designed as read-only, and changes in json invalidate the cache.

Example

model = load_json_model("some_model.json")
model.json # see the actual underlying JSON
reactions(model) # see the list of reactions
source
COBREXA.boundsMethod
bounds(model::JSONModel)

Get the bounds for reactions, assuming the information is stored in .lower_bound and .upper_bound.

source
COBREXA.objectiveMethod
objective(model::JSONModel)

Collect .objective_coefficient keys from model reactions.

source
COBREXA.reactionsMethod
reactions(model::JSONModel)

Extract reaction names (stored as .id) from JSON model.

source
COBREXA.stoichiometryMethod
stoichiometry(model::JSONModel)

Get the stoichiometry. Assuming the information is stored in reaction object under key .metabolites.

source
COBREXA.MATModelType
struct MATModel

Wrapper around the models loaded in dictionaries from the MATLAB representation.

source
Base.convertMethod
Base.convert(::Type{MATModel}, m::MetabolicModel)

Convert any metabolic model to MATModel.

source
COBREXA.balanceMethod
balance(m::MATModel)

Extracts balance from the MAT model, defaulting to zeroes if not present.

source
COBREXA.boundsMethod
bounds(m::MATModel)

Extracts bounds from the MAT file, saved under lb and ub.

source
COBREXA.coupling_boundsMethod
coupling_bounds(m::MATModel)

Extracts the coupling constraints. Currently, there are several accepted ways to store these in MATLAB models; this takes the constraints from vectors cl and cu.

source
COBREXA.genesMethod
genes(m::MATModel)

Extracts the possible gene list from genes key.

source
COBREXA.metabolitesMethod
metabolites(m::MATModel)::Vector{String}

Extracts metabolite names from mets key in the MAT file.

source
COBREXA.objectiveMethod
objective(m::MATModel)

Extracts the objective from the MAT model (defaults to zeroes).

source
COBREXA.reactionsMethod
reactions(m::MATModel)::Vector{String}

Extracts reaction names from rxns key in the MAT file.

source
COBREXA.MetaboliteType

Metabolite structure.

Fields

id :: String
name :: Maybe{String}
formula :: String
charge :: Int
compartment :: String
notes :: Dict{String, Vector{String}}
annotation :: Dict{String, Union{Vector{String}, String}}
source
COBREXA.ReactionType
Reaction(
    id::String,
    metabolites::Dict{String,Union{Int, Float64}},
    dir = :bidirectional;
    default_bound = _constants.default_reaction_bound,
)

Convenience constructor for Reaction. The reaction equation is specified using metabolites, which is a dictionary mapping metabolite ids to stoichiometric coefficients. The direcion of the reaction is set through dir which can take :bidirectional, :forward, and :reverse as values. Finally, the default_bound is the value taken to mean infinity in the context of constraint based models, often this is set to a very high flux value like 1000.

source
COBREXA.ReactionType
Reaction(
    id = "";
    name = nothing,
    metabolites = Dict{String,Float64}(),
    lb = -_constants.default_reaction_bound,
    ub = _constants.default_reaction_bound,
    grr = nothing,
    subsystem = nothing,
    notes = Notes(),
    annotations = Annotations(),
    objective_coefficient = 0.0,
)

A constructor for Reaction that only takes a reaction id and assigns default/uninformative values to all the fields that are not explicitely assigned.

source
COBREXA.ReactionType
mutable struct Reaction
    id::String
    name::Maybe{String}
    metabolites::Dict{String,Float64}
    lb::Float64
    ub::Float64
    grr::Maybe{GeneAssociation}
    subsystem::Maybe{String}
    notes::Notes
    annotations::Annotations
    objective_coefficient::Float64
end

A structure for representing a single reaction in a StandardModel.

source
COBREXA.SBMLModelType
struct SBMLModel

Thin wrapper around the model from SBML.jl library. Allows easy conversion from SBML to any other model format.

source
COBREXA.boundsMethod
bounds(model::SBMLModel)::Tuple{Vector{Float64},Vector{Float64}}

Get the lower and upper flux bounds of model SBMLModel. Throws DomainError in case if the SBML contains mismatching units.

source
COBREXA.SerializedType
mutable struct Serialized{M <: MetabolicModel}
    m::Maybe{M}
    filename::String
end

A meta-model that represents a model that is serialized on the disk. The internal model will be loaded on-demand by using any accessor, or by calling precache! directly.

source
COBREXA.precache!Method
precache!(model::Serialized{MetabolicModel})::Nothing

Load the Serialized model from disk in case it's not alreadly loaded.

source
COBREXA.StandardModelType
mutable struct StandardModel

StandardModel is used to store a constraint based metabolic model with meta-information. Meta-information is defined as annotation details, which include gene-reaction-rules, formulas, etc.

This model type seeks to keep as much meta-information as possible, as opposed to CoreModel and CoreModelCoupled, which keep the bare neccessities only. When merging models and keeping meta-information is important, use this as the model type. If meta-information is not important, use the more efficient core model types. See CoreModel and CoreModelCoupled for comparison.

In this model, reactions, metabolites, and genes are stored in ordered dictionaries indexed by each struct's id field. For example, model.reactions["rxn1_id"] returns a Reaction where the field id equals "rxn1_id". This makes adding and removing reactions efficient.

Note that the stoichiometric matrix (or any other core data, e.g. flux bounds) is not stored directly as in CoreModel. When this model type is used in analysis functions, these core data structures are built from scratch each time an analysis function is called. This can cause performance issues if you run many small analysis functions sequentially. Consider using the core model types if performance is critical.

See also: Reaction, Metabolite, Gene

Fields

id :: String
reactions :: OrderedDict{String, Reaction}
metabolites :: OrderedDict{String, Metabolite}
genes :: OrderedDict{String, Gene}

Example

model = load_model(StandardModel, "my_model.json")
keys(model.reactions)
source
Base.convertMethod

Base.convert(::Type{StandardModel}, model::MetabolicModel)

Convert any MetabolicModel into a StandardModel. Note, some data loss may occur since only the generic interface is used during the conversion process.

source
COBREXA.balanceMethod
balance(model::StandardModel)

Return the balance of the linear problem, i.e. b in Sv = 0 where S is the stoichiometric matrix and v is the flux vector.

source
COBREXA.boundsMethod
bounds(model::StandardModel)::Tuple{Vector{Float64},Vector{Float64}}

Return the lower and upper bounds, respectively, for reactions in model. Order matches that of the reaction ids returned in reactions().

source
COBREXA.gene_annotationsMethod
gene_annotations(model::StandardModel, id::String)::Annotations

Return the annotation associated with gene id in model. Return an empty Dict if not present.

source
COBREXA.gene_notesMethod
gene_notes(model::StandardModel, id::String)::Notes

Return the notes associated with gene id in model. Return an empty Dict if not present.

source
COBREXA.genesMethod
genes(model::StandardModel)

Return a vector of gene id strings in model.

source
COBREXA.lower_boundsMethod
lower_bounds(model::StandardModel)::Vector{Float64}

Return the lower bounds for all reactions in model in sparse format.

source
COBREXA.metabolite_annotationsMethod
metabolite_annotations(model::StandardModel, id::String)::Annotations

Return the annotation associated with metabolite id in model. Return an empty Dict if not present.

source
COBREXA.metabolite_chargeMethod
metabolite_charge(model::StandardModel, id::String)

Return the charge associated with metabolite id in model. Return nothing if not present.

source
COBREXA.metabolite_compartmentMethod
metabolite_compartment(model::StandardModel, id::String)

Return compartment associated with metabolite id in model. Return nothing if not present.

source
COBREXA.metabolite_formulaMethod
metabolite_formula(model::StandardModel, id::String)

Return the formula of reaction id in model. Return nothing if not present.

source
COBREXA.metabolite_notesMethod
metabolite_notes(model::StandardModel, id::String)::Notes

Return the notes associated with metabolite id in model. Return an empty Dict if not present.

source
COBREXA.metabolitesMethod
metabolites(model::StandardModel)

Return a vector of metabolite id strings contained in model. The order of metabolite strings returned here matches the order used to construct the stoichiometric matrix.

source
COBREXA.reaction_annotationsMethod
reaction_annotations(model::StandardModel, id::String)::Annotations

Return the annotation associated with reaction id in model. Return an empty Dict if not present.

source
COBREXA.reaction_gene_associationMethod
reaction_gene_association(model::StandardModel, id::String)

Return the gene reaction rule in string format for reaction with id in model. Return nothing if not available.

source
COBREXA.reaction_notesMethod
reaction_notes(model::StandardModel, id::String)::Notes

Return the notes associated with reaction id in model. Return an empty Dict if not present.

source
COBREXA.reaction_subsystemMethod
reaction_subsystem(id::String, model::StandardModel)

Return the subsystem associated with reaction id in model. Return nothing if not present.

source
COBREXA.reactionsMethod
reactions(model::StandardModel)

Return a vector of reaction id strings contained in model. The order of reaction ids returned here matches the order used to construct the stoichiometric matrix.

source
COBREXA.stoichiometryMethod
stoichiometry(model::StandardModel)

Return the stoichiometric matrix associated with model in sparse format.

source
COBREXA.upper_boundsMethod
upper_bounds(model::StandardModel)::Vector{Float64}

Return the upper bounds for all reactions in model in sparse format. Order matches that of the reaction ids returned in reactions().

source