Types
Base types
COBREXA.Maybe
— TypeMaybe{T} = Union{Nothing, T}
A nice name for "nullable" type.
COBREXA._default
— Method_default(d, x::Union{Nothing, T} where T) -> Any
Fold the Maybe{T}
down to T
by defaulting.
COBREXA._maybemap
— Method_maybemap(f, x) -> Any
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 of 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.ModelWrapper
— Typeabstract type ModelWrapper <: MetabolicModel end
A helper supertype of all "wrapper" types that contain precisely one other MetabolicModel
.
COBREXA.Notes
— TypeNotes = Dict{String,Vector{String}}
Free-form notes about something (e.g. a Gene
), categorized by "topic".
Model types and contents
COBREXA.CoreModel
— Typemutable 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ᵤ
Fields
S::SparseArrays.SparseMatrixCSC{Float64, Int64}
b::SparseArrays.SparseVector{Float64, Int64}
c::SparseArrays.SparseVector{Float64, Int64}
xl::Vector{Float64}
xu::Vector{Float64}
rxns::Vector{String}
mets::Vector{String}
grrs::Vector{Union{Nothing, Vector{Vector{String}}}}
Base.convert
— Methodconvert(_::Type{CoreModel}, m::MetabolicModel) -> MetabolicModel
Make a CoreModel
out of any compatible model type.
COBREXA.balance
— Methodbalance(a::CoreModel) -> SparseArrays.SparseVector{Float64, Int64}
CoreModel
target flux balance.
COBREXA.bounds
— Methodbounds(a::CoreModel) -> Tuple{Vector{Float64}, Vector{Float64}}
CoreModel
flux bounds.
COBREXA.genes
— Methodgenes(a::CoreModel) -> Vector{String}
Collect all genes contained in the CoreModel
. The call is expensive for large models, because the vector is not stored and instead gets rebuilt each time this function is called.
COBREXA.metabolites
— Methodmetabolites(a::CoreModel) -> Vector{String}
Metabolites in a CoreModel
.
COBREXA.objective
— Methodobjective(a::CoreModel) -> SparseArrays.SparseVector{Float64, Int64}
CoreModel
objective vector.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::CoreModel, ridx::Int64) -> Union{Nothing, Vector{Vector{String}}}
Retrieve the GeneAssociation
from CoreModel
by reaction index.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::CoreModel, rid::String) -> Union{Nothing, Vector{Vector{String}}}
Retrieve the GeneAssociation
from CoreModel
by reaction ID.
COBREXA.reaction_gene_association_vec
— Methodreaction_gene_association_vec(model::CoreModel) -> Vector{Union{Nothing, Vector{Vector{String}}}}
Retrieve a vector of all gene associations in a CoreModel
, in the same order as reactions(model)
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(m::CoreModel, ridx) -> Dict{String, Float64}
Return the stoichiometry of reaction at index ridx
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(m::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) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
CoreModel
stoichiometry matrix.
COBREXA.CoreCoupling
— Typemutable struct CoreCoupling{M} <: ModelWrapper
A matrix-based wrap that adds reaction coupling matrix to the inner model. A flux x
feasible in this model must satisfy:
cₗ ≤ C x ≤ cᵤ
Fields
lm::Any
C::SparseArrays.SparseMatrixCSC{Float64, Int64}
cl::Vector{Float64}
cu::Vector{Float64}
COBREXA.CoreModelCoupled
— Typeconst CoreModelCoupled = CoreCoupling{CoreModel}
A matrix-based linear model with additional coupling constraints in the form:
cₗ ≤ C x ≤ cᵤ
Internally, the model is implemented using CoreCoupling
that contains a single CoreModel
.
Base.convert
— Methodconvert(::Type{CoreCoupling{M}}, mm::MetabolicModel; clone_coupling) -> MetabolicModel
Make a CoreCoupling
out of any compatible model type.
COBREXA.coupling
— Methodcoupling(a::CoreCoupling) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
Coupling constraint matrix for a CoreCoupling
.
COBREXA.coupling_bounds
— Methodcoupling_bounds(a::CoreCoupling) -> Tuple{Vector{Float64}, Vector{Float64}}
Coupling bounds for a CoreCoupling
.
COBREXA.n_coupling_constraints
— Methodn_coupling_constraints(a::CoreCoupling) -> Int64
The number of coupling constraints in a CoreCoupling
.
COBREXA.reaction_gene_association_vec
— Methodreaction_gene_association_vec(model::CoreModelCoupled)
Evaluates reaction_gene_association_vec
on the model contained in CoreModelCoupled.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::CoreModelCoupled, ridx::Int)
Evaluates reaction_stoichiometry
on the model contained in CoreModelCoupled.
COBREXA.unwrap_model
— Methodunwrap_model(a::CoreCoupling) -> Any
Get the internal CoreModel
out of CoreCoupling
.
COBREXA.FluxSummary
— Typestruct FluxSummary
A struct used to store summary information about the solution of a constraint based analysis result.
Fields
biomass_fluxes::OrderedCollections.OrderedDict{String, Float64}
import_fluxes::OrderedCollections.OrderedDict{String, Float64}
export_fluxes::OrderedCollections.OrderedDict{String, Float64}
unbounded_fluxes::OrderedCollections.OrderedDict{String, Float64}
COBREXA.FluxSummary
— MethodFluxSummary() -> FluxSummary
A default empty constructor for FluxSummary
.
COBREXA.flux_summary
— Methodflux_summary(flux_result::Union{Nothing, Dict{String, Float64}}; exclude_exchanges, exchange_prefixes, biomass_strings, exclude_biomass, small_flux_bound, large_flux_bound, keep_unbounded) -> 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
— Typestruct FluxVariabilitySummary
Stores summary information about the result of a flux variability analysis.
Fields
biomass_fluxes::Dict{String, Vector{Union{Nothing, Float64}}}
exchange_fluxes::Dict{String, Vector{Union{Nothing, Float64}}}
COBREXA.FluxVariabilitySummary
— MethodFluxVariabilitySummary() -> FluxVariabilitySummary
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, exchange_prefixes, biomass_strings, exclude_biomass) -> 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() -> Gene
Gene(id; name, notes, annotations) -> Gene
A convenient constructor for a Gene
.
COBREXA.Gene
— Typemutable struct Gene
Fields
id::String
name::Union{Nothing, String}
notes::Dict{String, Vector{String}}
annotations::Dict{String, Vector{String}}
COBREXA.HDF5Model
— Typemutable struct HDF5Model <: MetabolicModel
A model that is stored in HDF5 format. The model data is never really pulled into memory, but instead mmap'ed as directly as possible into the Julia structures. This makes reading the HDF5Model
s extremely fast, at the same time the (uncached) HDF5Model
s can be sent around efficiently among distributed nodes just like Serialized
models, provided the nodes share a common storage.
All HDF5Models must have the backing disk storage. To create one, use save_h5_model
or save_model
with .h5
file extension. To create a temporary model that behaves like a model "in memory", save it to a temporary file. For related reasons, you can not use convert
models to HDF5Model
format, because the conversion would impliy having the model saved somewhere.
Fields
h5::Union{Nothing, HDF5.File}
filename::String
COBREXA.Isozyme
— Typemutable struct Isozyme
Information about isozyme composition and activity.
Fields
gene_product_count::Dict{String, Int64}
kcat_forward::Float64
kcat_reverse::Float64
COBREXA.JSONModel
— Typestruct JSONModel <: MetabolicModel
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
Fields
json::Dict{String, Any}
rxn_index::Dict{String, Int64}
rxns::Vector{Any}
met_index::Dict{String, Int64}
mets::Vector{Any}
gene_index::Dict{String, Int64}
genes::Vector{Any}
Base.convert
— Methodconvert(_::Type{JSONModel}, mm::MetabolicModel) -> MetabolicModel
Convert any MetabolicModel
to JSONModel
.
COBREXA.bounds
— Methodbounds(model::JSONModel) -> Tuple{Vector, Vector}
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) -> Dict{String, Vector{String}}
Gene annotations from the JSONModel
.
COBREXA.gene_name
— Methodgene_name(model::JSONModel, gid::String) -> Any
Return the name of gene with ID gid
.
COBREXA.gene_notes
— Methodgene_notes(model::JSONModel, gid::String) -> Dict{String, Vector{String}}
Gene notes from the JSONModel
.
COBREXA.genes
— Methodgenes(model::JSONModel) -> Vector
Extract gene names from a JSON model.
COBREXA.metabolite_annotations
— Methodmetabolite_annotations(model::JSONModel, mid::String) -> Dict{String, Vector{String}}
Metabolite annotations from the JSONModel
.
COBREXA.metabolite_charge
— Methodmetabolite_charge(model::JSONModel, mid::String) -> Any
Return the metabolite .charge
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(model::JSONModel, mid::String) -> Any
Return the metabolite .compartment
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::JSONModel, mid::String) -> Union{Nothing, Dict{String, Int64}}
Parse and return the metabolite .formula
COBREXA.metabolite_name
— Methodmetabolite_name(model::JSONModel, mid::String) -> Any
Return the name of metabolite with ID mid
.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::JSONModel, mid::String) -> Dict{String, Vector{String}}
Metabolite notes from the JSONModel
.
COBREXA.metabolites
— Methodmetabolites(model::JSONModel) -> Vector
Extract metabolite names (stored as .id
) from JSON model.
COBREXA.objective
— Methodobjective(model::JSONModel) -> SparseArrays.SparseVector{_A, Int64} where _A
Collect .objective_coefficient
keys from model reactions.
COBREXA.reaction_annotations
— Methodreaction_annotations(model::JSONModel, rid::String) -> Dict{String, Vector{String}}
Reaction annotations from the JSONModel
.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::JSONModel, rid::String) -> Union{Nothing, Vector{Vector{String}}}
Parses the .gene_reaction_rule
from reactions.
COBREXA.reaction_name
— Methodreaction_name(model::JSONModel, rid::String) -> Any
Return the name of reaction with ID rid
.
COBREXA.reaction_notes
— Methodreaction_notes(model::JSONModel, rid::String) -> Dict{String, Vector{String}}
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) -> Any
Parses the .subsystem
out from reactions.
COBREXA.reactions
— Methodreactions(model::JSONModel) -> Vector
Extract reaction names (stored as .id
) from JSON model.
COBREXA.stoichiometry
— Methodstoichiometry(model::JSONModel) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
Get the stoichiometry. Assuming the information is stored in reaction object under key .metabolites
.
COBREXA.MATModel
— Typestruct MATModel <: MetabolicModel
Wrapper around the models loaded in dictionaries from the MATLAB representation.
Fields
mat::Dict{String, Any}
Base.convert
— Methodconvert(_::Type{MATModel}, m::MetabolicModel) -> MetabolicModel
Convert any metabolic model to MATModel
.
COBREXA._mat_has_squashed_coupling
— Method_mat_has_squashed_coupling(mat) -> Any
Guesses whether C in the MAT file is stored in A=[S;C].
COBREXA.balance
— Methodbalance(m::MATModel) -> Any
Extracts balance from the MAT model, defaulting to zeroes if not present.
COBREXA.bounds
— Methodbounds(m::MATModel) -> Tuple{Any, Any}
Extracts bounds from the MAT file, saved under lb
and ub
.
COBREXA.coupling
— Methodcoupling(m::MATModel) -> Any
Extract coupling matrix stored, in C
key.
COBREXA.coupling_bounds
— Methodcoupling_bounds(m::MATModel) -> Tuple{Any, Any}
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) -> Any
Extracts the possible gene list from genes
key.
COBREXA.metabolite_charge
— Methodmetabolite_charge(m::MATModel, mid::String) -> Union{Nothing, Int64}
Extract metabolite charge from metCharge
or metCharges
.
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(m::MATModel, mid::String) -> Any
Extract metabolite compartment from metCompartment
or metCompartments
.
COBREXA.metabolite_formula
— Methodmetabolite_formula(m::MATModel, mid::String) -> Union{Nothing, Dict{String, Int64}}
Extract metabolite formula from key metFormula
or metFormulas
.
COBREXA.metabolite_name
— Methodmetabolite_name(m::MATModel, mid::String) -> Any
Extract metabolite name from metNames
.
COBREXA.metabolites
— Methodmetabolites(m::MATModel) -> Vector{String}
Extracts metabolite names from mets
key in the MAT file.
COBREXA.objective
— Methodobjective(m::MATModel) -> Any
Extracts the objective from the MAT model (defaults to zeroes).
COBREXA.reaction_gene_association
— Methodreaction_gene_association(m::MATModel, rid::String) -> Union{Nothing, Vector{Vector{String}}}
Extracts the associations from grRules
key, if present.
COBREXA.reaction_name
— Methodreaction_name(m::MATModel, rid::String) -> Any
Extract reaction name from rxnNames
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(m::MATModel, ridx) -> Dict{String, Float64}
Return the stoichiometry of reaction at index ridx
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(m::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) -> Any
Extract the stoichiometry matrix, stored under key S
.
COBREXA.balance
— Methodbalance(a::MetabolicModel) -> SparseArrays.SparseVector{Float64, Int64}
Get the sparse balance vector of a model.
COBREXA.bounds
— Methodbounds(a::MetabolicModel) -> Tuple{Vector{Float64}, Vector{Float64}}
Get the lower and upper solution bounds of a model.
COBREXA.coupling
— Methodcoupling(a::MetabolicModel) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
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.fluxes
— Methodfluxes(a::MetabolicModel) -> Vector{String}
In some models, the reactions
that correspond to the columns of stoichiometry
matrix do not fully represent the semantic contents of the model; for example, fluxes may be split into forward and reverse reactions, reactions catalyzed by distinct enzymes, etc. Together with reaction_flux
(and n_fluxes
) this specifies how the flux is decomposed into individual reactions.
By default (and in most models), fluxes and reactions perfectly correspond.
COBREXA.gene_annotations
— Methodgene_annotations(a::MetabolicModel, gene_id::String) -> Any
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) -> Any
Return the name of gene with ID gid
.
COBREXA.gene_notes
— Methodgene_notes(model::MetabolicModel, gene_id::String) -> Any
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) -> Any
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) -> Any
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) -> Any
Return the compartment of metabolite metabolite_id
in model
if it is assigned. If not, return nothing
.
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::MetabolicModel, metabolite_id::String) -> Any
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) -> Any
Return the name of metabolite with ID mid
.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::MetabolicModel, metabolite_id::String) -> Any
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. The vector precisely corresponds to the rows in stoichiometry
matrix.
As with reactions
s, some metabolites in models may be virtual, representing purely technical equality constraints.
COBREXA.n_coupling_constraints
— Methodn_coupling_constraints(a::MetabolicModel) -> Int64
Get the number of coupling constraints in a model.
COBREXA.n_genes
— Methodn_genes(a::MetabolicModel) -> Int64
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) -> Int64
Get the number of metabolites in a model.
COBREXA.n_reactions
— Methodn_reactions(a::MetabolicModel) -> Int64
Get the number of reactions in a model.
COBREXA.objective
— Methodobjective(a::MetabolicModel) -> SparseArrays.SparseVector{Float64, Int64}
Get the objective vector of the model. Analysis functions, such as flux_balance_analysis
, are supposed to maximize dot(objective, x)
where x
is a feasible solution of the model.
COBREXA.precache!
— Methodprecache!(a::MetabolicModel)
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) -> Any
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_flux
— Methodreaction_flux(a::MetabolicModel) -> Any
Retrieve a sparse matrix that describes the correspondence of a solution of the linear system to the fluxes (see fluxes
for rationale). Returns a sparse matrix of size (n_reactions(a), n_fluxes(a))
. For most models, this is an identity matrix.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(a::MetabolicModel, reaction_id::String) -> Any
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) -> Any
Return the name of reaction with ID rid
.
COBREXA.reaction_notes
— Methodreaction_notes(model::MetabolicModel, reaction_id::String) -> Any
Return the notes associated with reaction reaction_id
in model
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(m::MetabolicModel, rid::String) -> Any
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) -> Any
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. The vector precisely corresponds to the columns in stoichiometry
matrix.
For technical reasons, the "reactions" may sometimes not be true reactions but various virtual and helper pseudo-reactions that are used in the metabolic modeling, such as metabolite exchanges, separate forward and reverse reactions, supplies of enzymatic and genetic material and virtual cell volume, etc. To simplify the view of the model contents use reaction_flux
.
COBREXA.stoichiometry
— Methodstoichiometry(a::MetabolicModel) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
Get the sparse stoichiometry matrix of a model. A feasible solution x
of a model m
is defined as satisfying the equations:
stoichiometry(m) * x .== balance(m)
x .>= lbs
y .<= ubs
- `(lbs, ubs) == bounds(m)
COBREXA.Metabolite
— TypeMetabolite() -> Metabolite
Metabolite(id; name, formula, charge, compartment, notes, annotations) -> Metabolite
A constructor for Metabolite
s.
COBREXA.Metabolite
— Typemutable struct Metabolite
Fields
id::String
name::Union{Nothing, String}
formula::Union{Nothing, String}
charge::Union{Nothing, Int64}
compartment::Union{Nothing, String}
notes::Dict{String, Vector{String}}
annotations::Dict{String, Vector{String}}
COBREXA.balance
— Methodbalance(model::ModelWrapper)
Evaluates balance
on the model contained in ModelWrapper.
COBREXA.bounds
— Methodbounds(model::ModelWrapper)
Evaluates bounds
on the model contained in ModelWrapper.
COBREXA.coupling
— Methodcoupling(model::ModelWrapper)
Evaluates coupling
on the model contained in ModelWrapper.
COBREXA.coupling_bounds
— Methodcoupling_bounds(model::ModelWrapper)
Evaluates coupling_bounds
on the model contained in ModelWrapper.
COBREXA.fluxes
— Methodfluxes(model::ModelWrapper)
Evaluates fluxes
on the model contained in ModelWrapper.
COBREXA.gene_annotations
— Methodgene_annotations(model::ModelWrapper, gid::String)
Evaluates gene_annotations
on the model contained in ModelWrapper.
COBREXA.gene_notes
— Methodgene_notes(model::ModelWrapper, gid::String)
Evaluates gene_notes
on the model contained in ModelWrapper.
COBREXA.genes
— Methodgenes(model::ModelWrapper)
Evaluates genes
on the model contained in ModelWrapper.
COBREXA.metabolite_annotations
— Methodmetabolite_annotations(model::ModelWrapper, mid::String)
Evaluates metabolite_annotations
on the model contained in ModelWrapper.
COBREXA.metabolite_charge
— Methodmetabolite_charge(model::ModelWrapper, mid::String)
Evaluates metabolite_charge
on the model contained in ModelWrapper.
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(model::ModelWrapper, mid::String)
Evaluates metabolite_compartment
on the model contained in ModelWrapper.
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::ModelWrapper, mid::String)
Evaluates metabolite_formula
on the model contained in ModelWrapper.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::ModelWrapper, mid::String)
Evaluates metabolite_notes
on the model contained in ModelWrapper.
COBREXA.metabolites
— Methodmetabolites(model::ModelWrapper)
Evaluates metabolites
on the model contained in ModelWrapper.
COBREXA.n_coupling_constraints
— Methodn_coupling_constraints(model::ModelWrapper)
Evaluates n_coupling_constraints
on the model contained in ModelWrapper.
COBREXA.n_fluxes
— Methodn_fluxes(model::ModelWrapper)
Evaluates n_fluxes
on the model contained in ModelWrapper.
COBREXA.n_genes
— Methodn_genes(model::ModelWrapper)
Evaluates n_genes
on the model contained in ModelWrapper.
COBREXA.objective
— Methodobjective(model::ModelWrapper)
Evaluates objective
on the model contained in ModelWrapper.
COBREXA.precache!
— Methodprecache!(model::ModelWrapper)
Evaluates precache!
on the model contained in ModelWrapper.
COBREXA.reaction_annotations
— Methodreaction_annotations(model::ModelWrapper, rid::String)
Evaluates reaction_annotations
on the model contained in ModelWrapper.
COBREXA.reaction_flux
— Methodreaction_flux(model::ModelWrapper)
Evaluates reaction_flux
on the model contained in ModelWrapper.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::ModelWrapper, rid::String)
Evaluates reaction_gene_association
on the model contained in ModelWrapper.
COBREXA.reaction_notes
— Methodreaction_notes(model::ModelWrapper, rid::String)
Evaluates reaction_notes
on the model contained in ModelWrapper.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(model::ModelWrapper, rid::String)
Evaluates reaction_stoichiometry
on the model contained in ModelWrapper.
COBREXA.reaction_subsystem
— Methodreaction_subsystem(model::ModelWrapper, rid::String)
Evaluates reaction_subsystem
on the model contained in ModelWrapper.
COBREXA.reactions
— Methodreactions(model::ModelWrapper)
Evaluates reactions
on the model contained in ModelWrapper.
COBREXA.stoichiometry
— Methodstoichiometry(model::ModelWrapper)
Evaluates stoichiometry
on the model contained in ModelWrapper.
COBREXA.unwrap_model
— Methodunwrap_model(a::ModelWrapper) -> Any
A simple helper to pick the single w
COBREXA.Reaction
— TypeReaction(id::String, metabolites) -> Reaction
Reaction(id::String, metabolites, dir; default_bound) -> Reaction
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() -> Reaction
Reaction(id; name, metabolites, lb, ub, grr, subsystem, notes, annotations, objective_coefficient) -> Reaction
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
A structure for representing a single reaction in a StandardModel
.
Fields
id::String
name::Union{Nothing, String}
metabolites::Dict{String, Float64}
lb::Float64
ub::Float64
grr::Union{Nothing, Vector{Vector{String}}}
subsystem::Union{Nothing, String}
notes::Dict{String, Vector{String}}
annotations::Dict{String, Vector{String}}
objective_coefficient::Float64
COBREXA.ReactionStatus
— Typemutable struct ReactionStatus
Used for concise reporting of modeling results.
Fields
already_present::Bool
index::Int64
info::String
COBREXA.SBMLModel
— Typestruct SBMLModel <: MetabolicModel
Construct the SBML model and add the necessary cached indexes, possibly choosing an active objective.
COBREXA.SBMLModel
— Typestruct SBMLModel <: MetabolicModel
Thin wrapper around the model from SBML.jl library. Allows easy conversion from SBML to any other model format.
Fields
sbml::SBML.Model
reaction_ids::Vector{String}
reaction_idx::Dict{String, Int64}
metabolite_ids::Vector{String}
metabolite_idx::Dict{String, Int64}
gene_ids::Vector{String}
active_objective::String
Base.convert
— Methodconvert(_::Type{SBMLModel}, mm::MetabolicModel) -> MetabolicModel
Convert any metabolic model to SBMLModel
.
COBREXA.balance
— Methodbalance(model::SBMLModel) -> SparseArrays.SparseVector{Float64, Int64}
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_annotations
— Methodgene_annotations(model::SBMLModel, gid::String) -> Dict{String, Vector{String}}
Return the annotations of gene with ID gid
.
COBREXA.gene_name
— Methodgene_name(model::SBMLModel, gid::String) -> Union{Nothing, String}
Return the name of gene with ID gid
.
COBREXA.gene_notes
— Methodgene_notes(model::SBMLModel, gid::String) -> Dict{String, Vector{String}}
Return the notes about gene with ID gid
.
COBREXA.genes
— Methodgenes(model::SBMLModel) -> Vector{String}
Get genes of a SBMLModel
.
COBREXA.metabolite_annotations
— Methodmetabolite_annotations(model::SBMLModel, mid::String) -> Dict{String, Vector{String}}
Return the annotations of metabolite with ID mid
.
COBREXA.metabolite_charge
— Methodmetabolite_charge(model::SBMLModel, mid::String) -> Union{Nothing, Int64}
Get charge of a chosen metabolite from SBMLModel
.
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(model::SBMLModel, mid::String) -> String
Get the compartment of a chosen metabolite from SBMLModel
.
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::SBMLModel, mid::String) -> Union{Nothing, Dict{String, Int64}}
Get MetaboliteFormula
from a chosen metabolite from SBMLModel
.
COBREXA.metabolite_name
— Methodmetabolite_name(model::SBMLModel, mid::String) -> Union{Nothing, String}
Return the name of metabolite with ID mid
.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::SBMLModel, mid::String) -> Dict{String, Vector{String}}
Return the notes about metabolite with ID mid
.
COBREXA.metabolites
— Methodmetabolites(model::SBMLModel) -> Vector{String}
Get metabolites from a SBMLModel
.
COBREXA.n_genes
— Methodn_genes(model::SBMLModel) -> Int64
Get number of genes in SBMLModel
.
COBREXA.n_metabolites
— Methodn_metabolites(model::SBMLModel) -> Int64
Efficient counting of metabolites in SBMLModel
.
COBREXA.n_reactions
— Methodn_reactions(model::SBMLModel) -> Int64
Efficient counting of reactions in SBMLModel
.
COBREXA.objective
— Methodobjective(model::SBMLModel) -> SparseArrays.SparseVector{Float64, Int64}
Objective of the SBMLModel
.
COBREXA.reaction_annotations
— Methodreaction_annotations(model::SBMLModel, rid::String) -> Dict{String, Vector{String}}
Return the annotations of reaction with ID rid
.
COBREXA.reaction_gene_association
— Methodreaction_gene_association(model::SBMLModel, rid::String) -> Union{Nothing, Vector{Vector{String}}}
Retrieve the GeneAssociation
from SBMLModel
.
COBREXA.reaction_name
— Methodreaction_name(model::SBMLModel, rid::String) -> Union{Nothing, String}
Return the name of reaction with ID rid
.
COBREXA.reaction_notes
— Methodreaction_notes(model::SBMLModel, rid::String) -> Dict{String, Vector{String}}
Return the notes about reaction with ID rid
.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(m::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) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
Recreate the stoichiometry matrix from the SBMLModel
.
COBREXA.Serialized
— Typemutable struct Serialized{M} <: ModelWrapper
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.
Fields
m::Union{Nothing, M} where M
filename::String
COBREXA.precache!
— Methodprecache!(model::Serialized)
Load the Serialized
model from disk in case it's not alreadly loaded.
COBREXA.unwrap_model
— Methodunwrap_model(m::Serialized) -> Union{Nothing, M} where M
Unwrap the serialized model (precaching it transparently).
COBREXA.StandardModel
— Typemutable struct StandardModel <: MetabolicModel
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
Example
model = load_model(StandardModel, "my_model.json")
keys(model.reactions)
Fields
id::String
reactions::OrderedCollections.OrderedDict{String, Reaction}
metabolites::OrderedCollections.OrderedDict{String, Metabolite}
genes::OrderedCollections.OrderedDict{String, Gene}
Base.convert
— Methodconvert(_::Type{StandardModel}, model::MetabolicModel) -> 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) -> SparseArrays.SparseVector{Float64, Int64}
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) -> Dict{String, Vector{String}}
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) -> Union{Nothing, String}
Return the name of gene with ID id
.
COBREXA.gene_notes
— Methodgene_notes(model::StandardModel, id::String) -> Dict{String, Vector{String}}
Return the notes associated with gene id
in model
. Return an empty Dict if not present.
COBREXA.genes
— Methodgenes(model::StandardModel) -> Vector{String}
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) -> Dict{String, Vector{String}}
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) -> Union{Nothing, Int64}
Return the charge associated with metabolite id
in model
. Return nothing if not present.
COBREXA.metabolite_compartment
— Methodmetabolite_compartment(model::StandardModel, id::String) -> Union{Nothing, String}
Return compartment associated with metabolite id
in model
. Return nothing
if not present.
COBREXA.metabolite_formula
— Methodmetabolite_formula(model::StandardModel, id::String) -> Union{Nothing, Dict{String, Int64}}
Return the formula of reaction id
in model
. Return nothing
if not present.
COBREXA.metabolite_name
— Methodmetabolite_name(m::StandardModel, mid::String) -> Union{Nothing, String}
Return the name of metabolite with ID id
.
COBREXA.metabolite_notes
— Methodmetabolite_notes(model::StandardModel, id::String) -> Dict{String, Vector{String}}
Return the notes associated with metabolite id
in model
. Return an empty Dict if not present.
COBREXA.metabolites
— Methodmetabolites(model::StandardModel) -> Vector{String}
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) -> Int64
Return the number of genes in model
.
COBREXA.n_metabolites
— Methodn_metabolites(model::StandardModel) -> Int64
Return the number of metabolites in model
.
COBREXA.n_reactions
— Methodn_reactions(model::StandardModel) -> Int64
Return the number of reactions contained in model
.
COBREXA.objective
— Methodobjective(model::StandardModel) -> SparseArrays.SparseVector{Float64, Int64}
Return sparse objective vector for model
.
COBREXA.reaction_annotations
— Methodreaction_annotations(model::StandardModel, id::String) -> Dict{String, Vector{String}}
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) -> Union{Nothing, Vector{Vector{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) -> Union{Nothing, String}
Return the name of reaction with ID id
.
COBREXA.reaction_notes
— Methodreaction_notes(model::StandardModel, id::String) -> Dict{String, Vector{String}}
Return the notes associated with reaction id
in model
. Return an empty Dict if not present.
COBREXA.reaction_stoichiometry
— Methodreaction_stoichiometry(m::StandardModel, rid::String) -> Dict{String, Float64}
Return the stoichiometry of reaction with ID rid
.
COBREXA.reaction_subsystem
— Methodreaction_subsystem(model::StandardModel, id::String) -> Union{Nothing, String}
Return the subsystem associated with reaction id
in model
. Return nothing
if not present.
COBREXA.reactions
— Methodreactions(model::StandardModel) -> Vector{String}
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) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
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()
.
Model type wrappers
COBREXA.ExpressionLimitedModel
— Typestruct ExpressionLimitedModel <: ModelWrapper
ExpressionLimitedModel
follows the methodology of the E-Flux algorithm to constraint the flux through the reactions in order to simulate the limited expression of genes (and thus the limited availability of the gene products).
Use make_expression_limited_model
or with_expression_limits
to construct the models.
E-Flux algorithm is closer described by: Colijn, Caroline, Aaron Brandes, Jeremy Zucker, Desmond S. Lun, Brian Weiner, Maha R. Farhat, Tan-Yun Cheng, D. Branch Moody, Megan Murray, and James E. Galagan. "Interpreting expression data with metabolic flux models: predicting Mycobacterium tuberculosis mycolic acid production." PLoS computational biology 5, no. 8 (2009): e1000489.
Fields
relative_expression::Dict{String, Float64}
Relative gene expression w.r.t. to some chosen reference; the normalization and scale of the values should match the expectations of the
bounding_function
.
inner::MetabolicModel
The wrapped model.
bounding_function::Function
Function used to calculate the new reaction bounds from expression. In
make_expression_limited_model
this defaults toexpression_probabilistic_bounds
.
COBREXA.GeckoModel
— Typestruct GeckoModel <: ModelWrapper
A model with complex enzyme concentration and capacity bounds, as described in Sánchez, Benjamín J., et al. "Improving the phenotype predictions of a yeast genome-scale metabolic model by incorporating enzymatic constraints." Molecular systems biology 13.8 (2017): 935.
Use make_gecko_model
or with_gecko
to construct this kind of model.
The model wraps another "internal" model, and adds following modifications:
- enzymatic reactions with known enzyme information are split into multiple forward and reverse variants for each isozyme,
- reaction coupling is added to ensure the groups of isozyme reactions obey the global reaction flux bounds from the original model,
- gene concentrations specified by each reaction and its gene product stoichiometry, can constrained by the user to reflect measurements, such as from mass spectrometry,
- additional coupling is added to simulate total masses of different proteins grouped by type (e.g., membrane-bound and free-floating proteins), which can be again constrained by the user (this is slightly generalized from original GECKO algorithm, which only considers a single group of indiscernible proteins).
The structure contains fields columns
that describe the contents of the stoichiometry matrix columns, coupling_row_reaction
, coupling_row_gene_product
and coupling_row_mass_group
that describe correspondence of the coupling rows to original model and determine the coupling bounds (note: the coupling for gene product is actually added to stoichiometry, not in coupling
), and inner
, which is the original wrapped model. The objective
of the model includes also the extra columns for individual genes, as held by coupling_row_gene_product
.
Implementation exposes the split reactions (available as reactions(model)
), but retains the original "simple" reactions accessible by fluxes
. The related constraints are implemented using coupling
and coupling_bounds
.
Fields
objective::SparseArrays.SparseVector{Float64, Int64}
columns::Vector{COBREXA._gecko_reaction_column}
coupling_row_reaction::Vector{Int64}
coupling_row_gene_product::Vector{Tuple{Int64, Tuple{Float64, Float64}}}
coupling_row_mass_group::Vector{COBREXA._gecko_capacity}
inner::MetabolicModel
COBREXA._gecko_capacity
— Typestruct _gecko_capacity
A helper struct that contains the gene product capacity terms organized by the grouping type, e.g. metabolic or membrane groups etc.
Fields
group_id::String
gene_product_idxs::Vector{Int64}
gene_product_molar_masses::Vector{Float64}
group_upper_bound::Float64
COBREXA._gecko_reaction_column
— Typestruct _gecko_reaction_column
A helper type for describing the contents of GeckoModel
s.
Fields
reaction_idx::Int64
isozyme_idx::Int64
direction::Int64
reaction_coupling_row::Int64
lb::Float64
ub::Float64
gene_product_coupling::Vector{Tuple{Int64, Float64}}
COBREXA.balance
— Methodbalance(model::GeckoModel) -> Any
Return the balance of the reactions in the inner model, concatenated with a vector of zeros representing the enzyme balance of a GeckoModel
.
COBREXA.bounds
— Methodbounds(model::GeckoModel) -> Tuple{Vector{Float64}, Vector{Float64}}
Return variable bounds for GeckoModel
.
COBREXA.coupling
— Methodcoupling(model::GeckoModel) -> Any
Return the coupling of GeckoModel
. That combines the coupling of the wrapped model, coupling for split (arm) reactions, and the coupling for the total enzyme capacity.
COBREXA.coupling_bounds
— Methodcoupling_bounds(model::GeckoModel) -> Tuple{Any, Any}
The coupling bounds for GeckoModel
(refer to coupling
for details).
COBREXA.genes
— Methodgenes(model::GeckoModel) -> Any
Return the gene ids of genes that have enzymatic constraints associated with them.
COBREXA.metabolites
— Methodmetabolites(model::GeckoModel) -> Any
Return the ids of all metabolites, both real and pseudo, for a GeckoModel
.
COBREXA.n_coupling_constraints
— Methodn_coupling_constraints(model::GeckoModel) -> Any
Count the coupling constraints in GeckoModel
(refer to coupling
for details).
COBREXA.n_genes
— Methodn_genes(model::GeckoModel) -> Int64
Return the number of genes that have enzymatic constraints associated with them.
COBREXA.n_metabolites
— Methodn_metabolites(model::GeckoModel) -> Any
Return the number of metabolites, both real and pseudo, for a GeckoModel
.
COBREXA.n_reactions
— Methodn_reactions(model::GeckoModel) -> Int64
Returns the number of all irreversible reactions in model
as well as the number of gene products that take part in enzymatic reactions.
COBREXA.objective
— Methodobjective(model::GeckoModel) -> SparseArrays.SparseVector{Float64, Int64}
Return the objective of the GeckoModel
. Note, the objective is with respect to the internal variables, i.e. reactions(model)
, which are the unidirectional reactions and the genes involved in enzymatic reactions that have kinetic data.
COBREXA.reaction_flux
— Methodreaction_flux(model::GeckoModel) -> Any
Get the mapping of the reaction rates in GeckoModel
to the original fluxes in the wrapped model.
COBREXA.reactions
— Methodreactions(model::GeckoModel) -> Any
Returns the internal reactions in a GeckoModel
(these may be split to forward- and reverse-only parts with different isozyme indexes; reactions IDs are mangled accordingly with suffixes).
COBREXA.stoichiometry
— Methodstoichiometry(model::GeckoModel) -> Any
Return a stoichiometry of the GeckoModel
. The enzymatic reactions are split into unidirectional forward and reverse ones, each of which may have multiple variants per isozyme.
COBREXA.SMomentModel
— Typestruct SMomentModel <: ModelWrapper
An enzyme-capacity-constrained model using sMOMENT algorithm, as described by Bekiaris, Pavlos Stephanos, and Steffen Klamt, "Automatic construction of metabolic models with enzyme constraints" BMC bioinformatics, 2020.
Use make_smoment_model
or with_smoment
to construct the models.
The model is constructed as follows:
- stoichiometry of the original model is retained as much as possible, but enzymatic reations are split into forward and reverse parts (marked by a suffix like
...#forward
and...#reverse
), - coupling is added to simulate a virtual metabolite "enzyme capacity", which is consumed by all enzymatic reactions at a rate given by enzyme mass divided by the corresponding kcat,
- the total consumption of the enzyme capacity is constrained to a fixed maximum.
The SMomentModel
structure contains a worked-out representation of the optimization problem atop a wrapped MetabolicModel
, in particular the separation of certain reactions into unidirectional forward and reverse parts, an "enzyme capacity" required for each reaction, and the value of the maximum capacity constraint. Original coupling in the inner model is retained.
In the structure, the field columns
describes the correspondence of stoichiometry columns to the stoichiometry and data of the internal wrapped model, and total_enzyme_capacity
is the total bound on the enzyme capacity consumption as specified in sMOMENT algorithm.
This implementation allows easy access to fluxes from the split reactions (available in reactions(model)
), while the original "simple" reactions from the wrapped model are retained as fluxes
. All additional constraints are implemented using coupling
and coupling_bounds
.
Fields
columns::Vector{COBREXA._smoment_column}
total_enzyme_capacity::Float64
inner::MetabolicModel
COBREXA._smoment_column
— Typestruct _smoment_column
A helper type that describes the contents of SMomentModel
s.
Fields
reaction_idx::Int64
direction::Int64
lb::Float64
ub::Float64
capacity_required::Float64
COBREXA.bounds
— Methodbounds(model::SMomentModel) -> Tuple{Vector{Float64}, Vector{Float64}}
Return the variable bounds for SMomentModel
.
COBREXA.coupling
— Methodcoupling(model::SMomentModel) -> Any
Return the coupling of SMomentModel
. That combines the coupling of the wrapped model, coupling for split reactions, and the coupling for the total enzyme capacity.
COBREXA.coupling_bounds
— Methodcoupling_bounds(model::SMomentModel) -> Tuple{Any, Any}
The coupling bounds for SMomentModel
(refer to coupling
for details).
COBREXA.n_coupling_constraints
— Methodn_coupling_constraints(model::SMomentModel) -> Any
Count the coupling constraints in SMomentModel
(refer to coupling
for details).
COBREXA.n_reactions
— Methodn_reactions(model::SMomentModel) -> Int64
The number of reactions (including split ones) in SMomentModel
.
COBREXA.objective
— Methodobjective(model::SMomentModel) -> Any
Reconstruct an objective of the SMomentModel
.
COBREXA.reaction_flux
— Methodreaction_flux(model::SMomentModel) -> Any
Get the mapping of the reaction rates in SMomentModel
to the original fluxes in the wrapped model.
COBREXA.reactions
— Methodreactions(model::SMomentModel) -> Vector
Returns the internal reactions in a SMomentModel
(these may be split to forward- and reverse-only parts; reactions IDs are mangled accordingly with suffixes).
COBREXA.stoichiometry
— Methodstoichiometry(model::SMomentModel) -> Any
Return a stoichiometry of the SMomentModel
. The enzymatic reactions are split into unidirectional forward and reverse ones.