Types
Base types
COBREXA.Maybe
— TypeMaybe{T} = Union{Nothing, T}
A nice name for "nullable" type.
COBREXA._default
— Method_default(d::T, x::Maybe{T})::T where {T}
Fold the Maybe{T}
down to T
by defaulting.
COBREXA._maybemap
— Method_maybemap(f, x::Maybe)::Maybe
Apply a function to x
only if it is not nothing
.
COBREXA.Annotations
— TypeAnnotations = 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"])
COBREXA.GeneAssociation
— TypeGeneAssociation = 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.)
COBREXA.MetabolicModel
— Typeabstract 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.
COBREXA.MetaboliteFormula
— TypeMetaboliteFormula = Dict{String,Int}
Dictionary of atoms and their abundances in a molecule.
COBREXA.Notes
— TypeNotes = Dict{String,Vector{String}}
Free-form notes about something (e.g. a Gene
), categorized by "topic".
COBREXA.balance
— Methodbalance(a::MetabolicModel)::SparseVec
Get the sparse balance vector of a model (ie. the b
from S x = b
).
COBREXA.bounds
— Methodbounds(a::MetabolicModel)::Tuple{Vector{Float64},Vector{Float64}}
Get the lower and upper flux bounds of a model.
COBREXA.coupling
— Methodcoupling(a::MetabolicModel)::SparseMat
Get a matrix of coupling constraint definitions of a model. By default, there is no coupling in the models.
COBREXA.coupling_bounds
— Methodcoupling_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.
COBREXA.gene_annotations
— Methodgene_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"]
.
COBREXA.gene_name
— Methodgene_name(model::MetabolicModel, gid::String)
Return the name of gene with ID gid
.
COBREXA.gene_notes
— Methodgene_notes(model::MetabolicModel, gene_id::String)::Notes
Return the notes associated with the gene gene_id
in model
.
COBREXA.genes
— Methodgenes(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.
COBREXA.metabolite_annotations
— Methodmetabolite_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"]
.
COBREXA.metabolite_charge
— Methodmetabolite_charge(model::MetabolicModel, metabolite_id::String)::Maybe{Int}
Return the charge associated with metabolite metabolite_id
in model
. Returns nothing
if charge not present.
COBREXA.metabolite_compartment
— Methodmetabolite_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
.
COBREXA.metabolite_formula
— Methodmetabolite_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.
COBREXA.metabolite_name
— Methodmetabolite_name(model::MetabolicModel, mid::String)
Return the name of metabolite with ID mid
.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::MetabolicModel, metabolite_id::String)::Notes
Return the notes associated with metabolite reaction_id
in model
.
COBREXA.metabolites
— Methodmetabolites(a::MetabolicModel)::Vector{String}
Return a vector of metabolite identifiers in a model.
COBREXA.n_coupling_constraints
— Methodn_coupling_constraints(a::MetabolicModel)::Int
Get the number of coupling constraints in a model.
COBREXA.n_genes
— Methodn_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.
COBREXA.n_metabolites
— Methodn_metabolites(a::MetabolicModel)::Int
Get the number of metabolites in a model.
COBREXA.n_reactions
— Methodn_reactions(a::MetabolicModel)::Int
Get the number of reactions in a model.
COBREXA.objective
— Methodobjective(a::MetabolicModel)::SparseVec
Get the objective vector of a model.
COBREXA.precache!
— Methodprecache!(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.
COBREXA.reaction_annotations
— Methodreaction_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"]
.
COBREXA.reaction_gene_association
— Methodreaction_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 [[]]
.)
COBREXA.reaction_name
— Methodreaction_name(model::MetabolicModel, rid::String)
Return the name of reaction with ID rid
.
COBREXA.reaction_notes
— Methodreaction_notes(model::MetabolicModel, reaction_id::String)::Notes
Return the notes associated with reaction reaction_id
in model
.
COBREXA.reaction_stoichiometry
— Methodreaction_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.
COBREXA.reaction_subsystem
— Methodreaction_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
.
COBREXA.reactions
— Methodreactions(a::MetabolicModel)::Vector{String}
Return a vector of reaction identifiers in a model.
COBREXA.stoichiometry
— Methodstoichiometry(a::MetabolicModel)::SparseMat
Get the sparse stoichiometry matrix of a model.
Model types and contents
COBREXA.CoreModel
— Typestruct 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ᵤ
Base.convert
— MethodBase.convert(::Type{CoreModel}, m::M) where {M <: MetabolicModel}
Make a CoreModel
out of any compatible model type.
COBREXA.balance
— Methodbalance(a::CoreModel)::SparseVec
CoreModel
target flux balance.
COBREXA.bounds
— Methodbounds(a::CoreModel)::Tuple{Vector{Float64},Vector{Float64}}
CoreModel
flux bounds.
COBREXA.metabolites
— Methodmetabolites(a::CoreModel)::Vector{String}
Metabolites in a CoreModel
.
COBREXA.objective
— Methodobjective(a::CoreModel)::SparseVec
CoreModel
objective vector.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::CoreModel, ridx::Int)::Maybe{GeneAssociation}
Retrieve the GeneAssociation
from CoreModel
by reaction index.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::CoreModel, rid::String)::Maybe{GeneAssociation}
Retrieve the GeneAssociation
from CoreModel
by reaction ID.
COBREXA.reaction_gene_association_vec
— Methodreaction_gene_association_vec(model::CoreModel)::Vector{Maybe{GeneAssociation}}
Retrieve a vector of all gene associations in a CoreModel
, in the same order as reactions(model)
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::CoreModel, ridx)::Dict{String, Float64}
Return the stoichiometry of reaction at index ridx
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::CoreModel, rid::String)::Dict{String, Float64}
Return the stoichiometry of reaction with ID rid
.
COBREXA.reactions
— Methodreactions(a::CoreModel)::Vector{String}
Get the reactions in a CoreModel
.
COBREXA.stoichiometry
— Methodstoichiometry(a::CoreModel)::SparseMat
CoreModel
stoichiometry matrix.
COBREXA.CoreModelCoupled
— Typestruct CoreModelCoupled <: MetabolicModel
The linear model with additional coupling constraints in the form
cₗ ≤ C x ≤ cᵤ
Base.convert
— MethodBase.convert(::Type{CoreModelCoupled}, mm::MetabolicModel)
Make a CoreModelCoupled
out of any compatible model type.
COBREXA.balance
— Methodbalance(a::CoreModelCoupled)
Extract balance from CoreModelCoupled
(uses the internal CoreModel
).
COBREXA.bounds
— Methodbounds(a::CoreModelCoupled)
Extract bounds from CoreModelCoupled
(uses the internal CoreModel
).
COBREXA.coupling
— Methodcoupling(a::CoreModelCoupled)::SparseMat
Coupling constraint matrix for a CoreModelCoupled
.
COBREXA.coupling_bounds
— Methodcoupling_bounds(a::CoreModelCoupled)::Tuple{Vector{Float64},Vector{Float64}}
Coupling bounds for a CoreModelCoupled
.
COBREXA.metabolites
— Methodmetabolites(a::CoreModelCoupled)
Extract metabolites from CoreModelCoupled
(uses the internal CoreModel
).
COBREXA.n_coupling_constraints
— Methodn_coupling_constraints(a::CoreModelCoupled)::Int
The number of coupling constraints in a CoreModelCoupled
.
COBREXA.objective
— Methodobjective(a::CoreModelCoupled)
Extract objective from CoreModelCoupled
(uses the internal CoreModel
).
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::CoreModelCoupled, ridx::Int)::Maybe{GeneAssociation}
Retrieve the GeneAssociation
from CoreModelCoupled
by reaction index.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::CoreModelCoupled, rid::String)::Maybe{GeneAssociation}
Retrieve the GeneAssociation
from CoreModelCoupled
by reaction ID.
COBREXA.reaction_gene_association_vec
— Methodreaction_gene_association_vec(model::CoreModelCoupled)::Vector{Maybe{GeneAssociation}}
Retrieve a vector of gene associations in a CoreModelCoupled
, in the same order as reactions(model)
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::CoreModelCoupled, ridx)::Dict{String, Float64}
Return the stoichiometry of reaction at index ridx
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::CoreModelCoupled, rid::String)::Dict{String, Float64}
Return the stoichiometry of reaction with ID rid
.
COBREXA.reactions
— Methodreactions(a::CoreModelCoupled)
Extract reactions from CoreModelCoupled
(uses the internal CoreModel
).
COBREXA.stoichiometry
— Methodstoichiometry(a::CoreModelCoupled)
Extract stoichiometry from CoreModelCoupled
(uses the internal CoreModel
).
COBREXA.FluxSummary
— TypeFluxSummary
A struct used to store summary information about the solution of a constraint based analysis result.
COBREXA.FluxSummary
— MethodFluxSummary()
A default empty constructor for FluxSummary
.
COBREXA.flux_summary
— Methodflux_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
COBREXA.FluxVariabilitySummary
— TypeFluxVariabilitySummary
Stores summary information about the result of a flux variability analysis.
COBREXA.FluxVariabilitySummary
— MethodFluxVariabilitySummary()
A default empty constructor for FluxVariabilitySummary
.
COBREXA.flux_variability_summary
— Methodflux_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
... ... ...
COBREXA.Gene
— TypeGene struct.
Fields
id :: String
name :: Maybe{String}
notes :: Dict{String, Vector{String}}
annotation :: Dict{String, Union{Vector{String}, String}}
COBREXA.JSONModel
— Typestruct 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
Base.convert
— MethodBase.convert(::Type{JSONModel}, mm::MetabolicModel)
Convert any MetabolicModel
to JSONModel
.
COBREXA.bounds
— Methodbounds(model::JSONModel)
Get the bounds for reactions, assuming the information is stored in .lower_bound
and .upper_bound
.
COBREXA.gene_annotations
— Methodgene_annotations(model::JSONModel, gid::String)::Annotations
Gene annotations from the JSONModel
.
COBREXA.gene_name
— Methodgene_name(model::JSONModel, gid::String)
Return the name of gene with ID gid
.
COBREXA.gene_notes
— Methodgene_notes(model::JSONModel, gid::String)::Notes
Gene notes from the JSONModel
.
COBREXA.genes
— Methodgenes(model::JSONModel)
Extract gene names from a JSON model.
COBREXA.metabolite_annotations
— Methodmetabolite_annotations(model::JSONModel, mid::String)::Annotations
Metabolite annotations from the JSONModel
.
COBREXA.metabolite_charge
— Methodmetabolite_charge(model::JSONModel, mid::String)
Return the metabolite .charge
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(model::JSONModel, mid::String)
Return the metabolite .compartment
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::JSONModel, mid::String)
Parse and return the metabolite .formula
COBREXA.metabolite_name
— Methodmetabolite_name(model::JSONModel, mid::String)
Return the name of metabolite with ID mid
.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::JSONModel, mid::String)::Notes
Metabolite notes from the JSONModel
.
COBREXA.metabolites
— Methodmetabolites(model::JSONModel)
Extract metabolite names (stored as .id
) from JSON model.
COBREXA.objective
— Methodobjective(model::JSONModel)
Collect .objective_coefficient
keys from model reactions.
COBREXA.reaction_annotations
— Methodreaction_annotations(model::JSONModel, rid::String)::Annotations
Reaction annotations from the JSONModel
.
COBREXA.reaction_gene_association
— Methodreaction_gene_associaton(model::JSONModel, rid::String)
Parses the .gene_reaction_rule
from reactions.
COBREXA.reaction_name
— Methodreaction_name(model::JSONModel, rid::String)
Return the name of reaction with ID rid
.
COBREXA.reaction_notes
— Methodreaction_notes(model::JSONModel, rid::String)::Notes
Reaction notes from the JSONModel
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::JSONModel, rid::String)::Dict{String, Float64}
Return the stoichiometry of reaction with ID rid
.
COBREXA.reaction_subsystem
— Methodreaction_subsystem(model::JSONModel, rid::String)
Parses the .subsystem
out from reactions.
COBREXA.reactions
— Methodreactions(model::JSONModel)
Extract reaction names (stored as .id
) from JSON model.
COBREXA.stoichiometry
— Methodstoichiometry(model::JSONModel)
Get the stoichiometry. Assuming the information is stored in reaction object under key .metabolites
.
COBREXA.MATModel
— Typestruct MATModel
Wrapper around the models loaded in dictionaries from the MATLAB representation.
Base.convert
— MethodBase.convert(::Type{MATModel}, m::MetabolicModel)
Convert any metabolic model to MATModel
.
COBREXA._mat_has_squashed_coupling
— Method_mat_has_squashed_coupling(mat)
Guesses whether C in the MAT file is stored in A=[S;C].
COBREXA.balance
— Methodbalance(m::MATModel)
Extracts balance from the MAT model, defaulting to zeroes if not present.
COBREXA.bounds
— Methodbounds(m::MATModel)
Extracts bounds from the MAT file, saved under lb
and ub
.
COBREXA.coupling
— Methodcoupling(m::MATModel)
Extract coupling matrix stored, in C
key.
COBREXA.coupling_bounds
— Methodcoupling_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
.
COBREXA.genes
— Methodgenes(m::MATModel)
Extracts the possible gene list from genes
key.
COBREXA.metabolite_charge
— Methodmetabolite_charge(m::MATModel, mid::String)
Extract metabolite charge from metCharge
or metCharges
.
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(m::MATModel, mid::String)
Extract metabolite compartment from metCompartment
or metCompartments
.
COBREXA.metabolite_formula
— Methodmetabolite_formula(m::MATModel, mid::String)
Extract metabolite formula from key metFormula
or metFormulas
.
COBREXA.metabolites
— Methodmetabolites(m::MATModel)::Vector{String}
Extracts metabolite names from mets
key in the MAT file.
COBREXA.objective
— Methodobjective(m::MATModel)
Extracts the objective from the MAT model (defaults to zeroes).
COBREXA.reaction_gene_association
— Methodreaction_gene_association(m::MATModel, rid::String)
Extracts the associations from grRules
key, if present.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::MATModel, ridx)::Dict{String, Float64}
Return the stoichiometry of reaction at index ridx
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::MATModel, rid::String)::Dict{String, Float64}
Return the stoichiometry of reaction with ID rid
.
COBREXA.reactions
— Methodreactions(m::MATModel)::Vector{String}
Extracts reaction names from rxns
key in the MAT file.
COBREXA.stoichiometry
— Methodstoichiometry(m::MATModel)
Extract the stoichiometry matrix, stored under key S
.
COBREXA.Metabolite
— TypeMetabolite 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}}
COBREXA.Reaction
— TypeReaction(
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.
COBREXA.Reaction
— TypeReaction(
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.
COBREXA.Reaction
— Typemutable 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
.
COBREXA.ReactionStatus
— TypeUsed for concise reporting of modeling results.
COBREXA.SBMLModel
— Typestruct SBMLModel
Thin wrapper around the model from SBML.jl library. Allows easy conversion from SBML to any other model format.
Base.convert
— MethodBase.convert(::Type{SBMLModel}, mm::MetabolicModel)
Convert any metabolic model to SBMLModel
.
COBREXA.balance
— Methodbalance(model::SBMLModel)::SparseVec
Balance vector of a SBMLModel
. This is always zero.
COBREXA.bounds
— Methodbounds(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.
COBREXA.gene_name
— Methodgene_name(model::SBMLModel, gid::String)
Return the name of gene with ID gid
.
COBREXA.genes
— Methodgenes(model::SBMLModel)::Vector{String}
Get genes of a SBMLModel
.
COBREXA.metabolite_charge
— Methodmetabolite_charge(model::SBMLModel, mid::String)::Maybe{Int}
Get charge of a chosen metabolite from SBMLModel
.
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::SBMLModel, mid::String)::Maybe{MetaboliteFormula}
Get MetaboliteFormula
from a chosen metabolite from SBMLModel
.
COBREXA.metabolite_name
— Methodmetabolite_name(model::SBMLModel, mid::String)
Return the name of metabolite with ID mid
.
COBREXA.metabolites
— Methodmetabolites(model::SBMLModel)::Vector{String}
Get metabolites from a SBMLModel
.
COBREXA.n_genes
— Methodn_genes(model::SBMLModel)::Int
Get number of genes in SBMLModel
.
COBREXA.n_metabolites
— Methodn_metabolites(model::SBMLModel)::Int
Efficient counting of metabolites in SBMLModel
.
COBREXA.n_reactions
— Methodn_reactions(model::SBMLModel)::Int
Efficient counting of reactions in SBMLModel
.
COBREXA.objective
— Methodobjective(model::SBMLModel)::SparseVec
Objective of the SBMLModel
.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::SBMLModel, rid::String)::Maybe{GeneAssociation}
Retrieve the GeneAssociation
from SBMLModel
.
COBREXA.reaction_name
— Methodreaction_name(model::SBMLModel, rid::String)
Return the name of reaction with ID rid
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::SBMLModel, rid::String)::Dict{String, Float64}
Return the stoichiometry of reaction with ID rid
.
COBREXA.reactions
— Methodreactions(model::SBMLModel)::Vector{String}
Get reactions from a SBMLModel
.
COBREXA.stoichiometry
— Methodstoichiometry(model::SBMLModel)::SparseMat
Recreate the stoichiometry matrix from the SBMLModel
.
COBREXA.Serialized
— Typemutable 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.
COBREXA.precache!
— Methodprecache!(model::Serialized{MetabolicModel})::Nothing
Load the Serialized
model from disk in case it's not alreadly loaded.
COBREXA.StandardModel
— Typemutable 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)
Base.convert
— MethodBase.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.
COBREXA.balance
— Methodbalance(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.
COBREXA.bounds
— Methodbounds(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()
.
COBREXA.gene_annotations
— Methodgene_annotations(model::StandardModel, id::String)::Annotations
Return the annotation associated with gene id
in model
. Return an empty Dict if not present.
COBREXA.gene_name
— Methodgene_name(m::StandardModel, gid::String)
Return the name of gene with ID id
.
COBREXA.gene_notes
— Methodgene_notes(model::StandardModel, id::String)::Notes
Return the notes associated with gene id
in model
. Return an empty Dict if not present.
COBREXA.genes
— Methodgenes(model::StandardModel)
Return a vector of gene id strings in model
.
COBREXA.lower_bounds
— Methodlower_bounds(model::StandardModel)::Vector{Float64}
Return the lower bounds for all reactions in model
in sparse format.
COBREXA.metabolite_annotations
— Methodmetabolite_annotations(model::StandardModel, id::String)::Annotations
Return the annotation associated with metabolite id
in model
. Return an empty Dict if not present.
COBREXA.metabolite_charge
— Methodmetabolite_charge(model::StandardModel, id::String)
Return the charge associated with metabolite id
in model
. Return nothing if not present.
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(model::StandardModel, id::String)
Return compartment associated with metabolite id
in model
. Return nothing
if not present.
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::StandardModel, id::String)
Return the formula of reaction id
in model
. Return nothing
if not present.
COBREXA.metabolite_name
— Methodmetabolite_name(m::StandardModel, mid::String)
Return the name of metabolite with ID id
.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::StandardModel, id::String)::Notes
Return the notes associated with metabolite id
in model
. Return an empty Dict if not present.
COBREXA.metabolites
— Methodmetabolites(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.
COBREXA.n_genes
— Methodn_genes(model::StandardModel)
Return the number of genes in model
.
COBREXA.n_metabolites
— Methodn_metabolites(model::StandardModel)
Return the number of metabolites in model
.
COBREXA.n_reactions
— Methodn_reactions(model::StandardModel)
Return the number of reactions contained in model
.
COBREXA.objective
— Methodobjective(model::StandardModel)
Return sparse objective vector for model
.
COBREXA.reaction_annotations
— Methodreaction_annotations(model::StandardModel, id::String)::Annotations
Return the annotation associated with reaction id
in model
. Return an empty Dict if not present.
COBREXA.reaction_gene_association
— Methodreaction_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.
COBREXA.reaction_name
— Methodreaction_name(m::StandardModel, rid::String)
Return the name of reaction with ID id
.
COBREXA.reaction_notes
— Methodreaction_notes(model::StandardModel, id::String)::Notes
Return the notes associated with reaction id
in model
. Return an empty Dict if not present.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::StandardModel, rid::String)::Dict{String, Float64}
Return the stoichiometry of reaction with ID rid
.
COBREXA.reaction_subsystem
— Methodreaction_subsystem(id::String, model::StandardModel)
Return the subsystem associated with reaction id
in model
. Return nothing
if not present.
COBREXA.reactions
— Methodreactions(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.
COBREXA.stoichiometry
— Methodstoichiometry(model::StandardModel)
Return the stoichiometric matrix associated with model
in sparse format.
COBREXA.upper_bounds
— Methodupper_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()
.