Functions

Base Types

COBREXA.MaybeType
Maybe{T} = Union{Nothing, T}

A nice name for "nullable" type.

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

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

COBREXA._maybemapMethod
_maybemap(f, x::Maybe)::Maybe

Apply a function to x only if it is not nothing.

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"])
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.)

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.

COBREXA.MetaboliteFormulaType
MetaboliteFormula = Dict{String,Int}

Dictionary of atoms and their abundances in a molecule.

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

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

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

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

COBREXA.boundsMethod
bounds(a::MetabolicModel)::Tuple{SparseVec,SparseVec}

Get the lower and upper flux bounds of a model.

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.

COBREXA.coupling_boundsMethod
coupling_bounds(a::MetabolicModel)::Tuple{SparseVec,SparseVec}

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_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"].

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

Return the notes associated with the gene gene_id in model.

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.

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"].

COBREXA.metabolite_chargeMethod

metabolitecharge(model::MetabolicModel, metaboliteid::String)::Maybe{Int}

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

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.

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.

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

Return the notes associated with metabolite reaction_id in model.

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

Return a vector of metabolite identifiers in a model.

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.

COBREXA.n_metabolitesMethod
n_metabolites(a::MetabolicModel)::Int

Get the number of metabolites in a model.

COBREXA.n_reactionsMethod
n_reactions(a::MetabolicModel)::Int

Get the number of reactions in a model.

COBREXA.objectiveMethod
objective(a::MetabolicModel)::SparseVec

Get the objective vector of a model.

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.

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"].

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 [[]].)

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

Return the notes associated with reaction reaction_id in model.

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.

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.

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

Return a vector of reaction identifiers in a model.

COBREXA.stoichiometryMethod
stoichiometry(a::MetabolicModel)::SparseMat

Get the sparse stoichiometry matrix of a model.

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ᵤ
Base.convertMethod
Base.convert(::Type{CoreModel}, m::M) where {M <: MetabolicModel}

Make a CoreModel out of any compatible model type.

COBREXA.balanceMethod
balance(a::CoreModel)::SparseVec

CoreModel target flux balance.

COBREXA.boundsMethod
bounds(a::CoreModel)::Tuple{SparseVec,SparseVec}

CoreModel flux bounds.

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

Metabolites in a CoreModel.

COBREXA.objectiveMethod
objective(a::CoreModel)::SparseVec

CoreModel objective vector.

COBREXA.reaction_stoichiometryMethod
reaction_stoichiometry(model::CoreModel, ridx)::Dict{String, Float64}

Return the stoichiometry of reaction at index ridx.

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

Return the stoichiometry of reaction with ID rid.

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

Get the reactions in a CoreModel.

COBREXA.CoreModelCoupledType
struct CoreModelCoupled <: MetabolicModel

The linear model with additional coupling constraints in the form

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

Make a CoreModelCoupled out of any compatible model type.

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

Coupling constraint matrix for a CoreModelCoupled.

COBREXA.coupling_boundsMethod
coupling_bounds(a::CoreModelCoupled)::Tuple{SparseVec,SparseVec}

Coupling bounds for a CoreModelCoupled.

COBREXA.reaction_stoichiometryMethod
reaction_stoichiometry(model::CoreModelCoupled, ridx)::Dict{String, Float64}

Return the stoichiometry of reaction at index ridx.

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

Return the stoichiometry of reaction with ID rid.

COBREXA.FluxSummaryType
FluxSummary

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

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_balance_analysis_dict(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.FluxVariabilitySummaryType
FluxVariabilitySummary

A struct used to store summary information about the solution of a flux variability analysis result.

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 fluxvariabilityanalysisdict. The simplified summary representation is useful for pretty-printing and easily showing the most important results. Internally this function uses [`lookslikebiomassreaction](@ref) and [lookslikeexchange_reaction`](@ref). The corresponding keyword arguments passed to these functions. Use this if your model has non-standard ids for 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.GeneType

Gene struct.

Fields

id :: String
notes :: Dict{String, Vector{String}}
annotation :: Dict{String, Union{Vector{String}, String}}
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
COBREXA.boundsMethod
bounds(model::JSONModel)

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

COBREXA.genesMethod
genes(model::JSONModel)

Extract gene names from a JSON model.

COBREXA.metabolitesMethod
metabolites(model::JSONModel)

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

COBREXA.objectiveMethod
objective(model::JSONModel)

Collect .objective_coefficient keys from model reactions.

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

Return the stoichiometry of reaction with ID rid.

COBREXA.reactionsMethod
reactions(model::JSONModel)

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

COBREXA.stoichiometryMethod
stoichiometry(model::JSONModel)

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

COBREXA.MATModelType
struct MATModel

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

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

Convert any metabolic model to MATModel.

COBREXA.balanceMethod
balance(m::MATModel)

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

COBREXA.boundsMethod
bounds(m::MATModel)

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

COBREXA.couplingMethod
coupling(m::MATModel)

Extract coupling matrix stored, in C key.

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.

COBREXA.genesMethod
genes(m::MATModel)

Extracts the possible gene list from genes key.

COBREXA.metabolite_chargeMethod
metabolite_charge(m::MATModel, mid::String)

Extract metabolite charge from metCharge or metCharges.

COBREXA.metabolite_compartmentMethod
metabolite_compartment(m::MATModel, mid::String)

Extract metabolite compartment from metCompartment or metCompartments.

COBREXA.metabolite_formulaMethod
metabolite_formula(m::MATModel, mid::String)

Extract metabolite formula from key metFormula or metFormulas.

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

Extracts metabolite names from mets key in the MAT file.

COBREXA.objectiveMethod
objective(m::MATModel)

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

COBREXA.reaction_stoichiometryMethod
reaction_stoichiometry(model::MATModel, ridx)::Dict{String, Float64}

Return the stoichiometry of reaction at index ridx.

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

Return the stoichiometry of reaction with ID rid.

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

Extracts reaction names from rxns key in the MAT file.

COBREXA.stoichiometryMethod
stoichiometry(m::MATModel)

Extract the stoichiometry matrix, stored under key S.

COBREXA.MetaboliteType

Metabolite structure.

Fields

id :: String
formula :: String
charge :: Int
compartment :: String
notes :: Dict{String, Vector{String}}
annotation :: Dict{String, Union{Vector{String}, String}}
COBREXA.ReactionType
mutable struct Reaction
    id::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.ReactionType
Reaction(
    id = "";
    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.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.

COBREXA.SBMLModelType
struct SBMLModel

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

Base.convertMethod
Base.convert(::Type{SBMLModel}, mm::MetabolicModel)

Convert any metabolic model to SBMLModel.

COBREXA.boundsMethod
bounds(model::SBMLModel)::Tuple{SparseVec,SparseVec}

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

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

Return the stoichiometry of reaction with ID rid.

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.

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

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

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

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.

COBREXA.boundsMethod
bounds(model::StandardModel)

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

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.

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.

COBREXA.genesMethod
genes(model::StandardModel)

Return a vector of gene id strings in model.

COBREXA.lower_boundsMethod
lower_bounds(model::StandardModel)

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

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.

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

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

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

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

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

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

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.

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.

COBREXA.n_genesMethod
n_genes(model::StandardModel)

Return the number of genes in model.

COBREXA.n_metabolitesMethod

n_metabolites(model::StandardModel)

Return the number of metabolites in model.

COBREXA.n_reactionsMethod
n_reactions(model::StandardModel)

Return the number of reactions contained in model.

COBREXA.objectiveMethod
objective(model::StandardModel)

Return sparse objective vector for model.

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.

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.

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.

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

Return the stoichiometry of reaction with ID rid.

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

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

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.

COBREXA.stoichiometryMethod
stoichiometry(model::StandardModel)

Return the stoichiometric matrix associated with model in sparse format.

COBREXA.upper_boundsMethod
upper_bounds(model::StandardModel)

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

Base functions

COBREXA._constantsConstant

A named tuple that contains the magic values that are used globally for whatever purposes.

COBREXA.is_solvedMethod
is_solved(optmodel)

Return true if optmodel solved successfully (solution is optimal or locally optimal). Return false if any other termination status is reached. Termination status is defined in the documentation of JuMP.

COBREXA.make_optimization_modelMethod
make_optimization_model(
    model::MetabolicModel,
    optimizer;
    sense = MOI.MAX_SENSE,
)

Convert MetabolicModels to a JuMP model, place objectives and the equality constraint.

COBREXA.optimize_objectiveMethod
optimize_objective(optmodel)::Union{Float64,Nothing}

Shortcut for running JuMP optimize! on a model and returning the objective value, if solved.

COBREXA.set_optmodel_bound!Method
set_optmodel_bound!(vidx, opt_model;
    ub::Maybe{Real} = nothing,
    lb::Maybe{Real} = nothing,
)

Helper function to set the bounds of a variable in the model. Internally calls set_normalized_rhs from JuMP. If the bounds are set to nothing, they will not be changed.

File I/O and serialization

COBREXA.load_modelMethod
load_model(file_name::String)::MetabolicModel

Generic function for loading models that chooses a specific loader function from the file_name extension, or throws an error.

Currently, these model types are supported:

COBREXA.load_modelMethod
load_model(type::Type{T}, file_name::String)::T where T

Helper function tht loads the model using load_model and return it converted to type.

Example:

load_model(CoreModel, "mySBMLModel.xml")
COBREXA.save_modelMethod
save_model(model::MetabolicModel, file_name::String)

Generic function for saving models that chooses a specific writer function from the file_name extension, or throws an error.

Currently, these model types are supported:

COBREXA.load_json_modelMethod
load_json_model(filename::String)::JSONModel

Load and return a JSON-formatted model that is stored in file_name.

COBREXA.save_json_modelMethod
save_json_model(model::MetabolicModel, file_name::String)

Save a JSONModel in model to a JSON file file_name.

In case the model is not JSONModel, it will be converted automatically.

COBREXA.load_mat_modelMethod
load_mat_model(file_name::String)

Load and return a MATLAB file file_name that contains a COBRA-compatible model.

COBREXA.save_mat_modelMethod
save_mat_model(model::MetabolicModel, file_name::String; model_name::String="model")

Save a MATModel in model to a MATLAB file file_name in a format compatible with other MATLAB-based COBRA software.

In case the model is not MATModel, it will be converted automatically.

model_name is the identifier name for the whole model written to the MATLAB file; defaults to just "model".

COBREXA.load_sbml_modelMethod
load_sbml_model(file_name::String)::SBMLModel

Load and return a SBML XML model in file_name.

Pretty printing

Base.showMethod

Pretty printing of everything metabolic-modelish.

Base.showMethod
Base.show(io::IO, ::MIME"text/plain", m::Serialized{M}) where {M}

Show the Serialized model without unnecessarily loading it.

COBREXA._pretty_print_keyvalsMethod
_pretty_print_keyvals(
    io,
    def::String,
    payload::Dict
)

Specialization of _pretty_print_keyvals for dictionaries.

COBREXA._pretty_print_keyvalsMethod
_pretty_print_keyvals(
    io,
    def::String,
    payload::String
)

Specialization of _pretty_print_keyvals for plain strings.

Model reconstruction

COBREXA.add_reaction!Method
add_reaction!(model::CoreModel, rxn::Reaction)

Add rxn to model. The model must already contain the metabolites used by rxn in the model.

COBREXA.add_reactions!Method
add_reactions!(model::CoreModel, rxns::Vector{Reaction})

Add rxns to model efficiently. The model must already contain the metabolites used by rxns in the model.

COBREXA.add_reactionsMethod
add_reactions(m1::CoreModel, m2::CoreModel; check_consistency = false)

Add all reactions from m2 to m1.

COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModel,
    s::V1,
    b::V2,
    c::AbstractFloat,
    xl::AbstractFloat,
    xu::AbstractFloat,
    rxn::String,
    mets::K;
    check_consistency = false,
) where {V1<:VecType,V2<:VecType,K<:StringVecType}
COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModel,
    Sp::M,
    b::V,
    c::V,
    xl::V,
    xu::V,
    rxns::K,
    mets::K;
    check_consistency = false,
) where {M<:MatType,V<:VecType,K<:StringVecType}
COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModel,
    s::V1,
    b::V2,
    c::AbstractFloat,
    xl::AbstractFloat,
    xu::AbstractFloat;
    check_consistency = false,
) where {V1<:VecType,V2<:VecType}

Add reaction(s) to a CoreModel model m.

COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModel,
    Sp::M,
    b::V,
    c::V,
    xl::V,
    xu::V;
    check_consistency = false,
) where {M<:MatType,V<:VecType}
COBREXA.verify_consistencyMethod
verify_consistency(
    m::CoreModel,
    Sp::M,
    b::V,
    c::V,
    xl::V,
    xu::V,
    names::K,
    mets::K,
    new_reactions,
    new_metabolites,
) where {M<:MatType,V<:VecType,K<:StringVecType}

Check the consistency of given reactions with existing reactions in m.

TODO: work in progress, doesn't return consistency status.

COBREXA.add_coupling_constraints!Method
add_coupling_constraints!(
    m::CoreModelCoupled,
    c::VecType,
    cl::AbstractFloat,
    cu::AbstractFloat,
)

Overload for adding a single coupling constraint.

COBREXA.add_coupling_constraints!Method
add_coupling_constraints!(
    m::CoreModelCoupled,
    C::MatType,
    cl::V,
    cu::V,
) where {V<:VecType}

In-place add a single coupling constraint in form

    cₗ ≤ C x ≤ cᵤ
COBREXA.add_reactionsMethod
add_reactions(m1::CoreModelCoupled, m2::CoreModel; check_consistency = false)

Add all reactions from m2 to m1.

COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModelCoupled,
    s::V1,
    b::V2,
    c::AbstractFloat,
    xl::AbstractFloat,
    xu::AbstractFloat,
    rxn::String,
    mets::K;
    check_consistency = false,
) where {V1<:VecType,V2<:VecType,K<:StringVecType}
COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModelCoupled,
    Sp::M,
    b::V,
    c::V,
    xl::V,
    xu::V,
    rxns::K,
    mets::K;
    check_consistency = false,
) where {M<:MatType,V<:VecType,K<:StringVecType}
COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModelCoupled,
    s::V1,
    b::V2,
    c::AbstractFloat,
    xl::AbstractFloat,
    xu::AbstractFloat;
    check_consistency = false,
) where {V1<:VecType,V2<:VecType}

Add reaction(s) to a CoreModelCoupled model m.

COBREXA.add_reactionsMethod
add_reactions(
    m::CoreModelCoupled,
    Sp::M,
    b::V,
    c::V,
    xl::V,
    xu::V;
    check_consistency = false,
) where {M<:MatType,V<:VecType}
COBREXA.change_coupling_bounds!Method
change_coupling_bounds!(
    model::CoreModelCoupled,
    constraints::Vector{Int};
    cl::V = Float64[],
    cu::V = Float64[],
) where {V<:VecType}

Change the lower and/or upper bounds (cl and cu) for the given list of coupling constraints.

COBREXA.:←Method
←(
    substrates::Union{
        Nothing,
        Metabolite,
        MetaboliteWithCoefficient,
        Vector{MetaboliteWithCoefficient},
    },
    products::Union{
        Nothing,
        Metabolite,
        MetaboliteWithCoefficient,
        Vector{MetaboliteWithCoefficient}
    },
)

Make a reverse-only Reaction from substrates and products.

COBREXA.:→Method
→(
    substrates::Union{
        Nothing,
        Metabolite,
        MetaboliteWithCoefficient,
        Vector{MetaboliteWithCoefficient},
    },
    products::Union{
        Nothing,
        Metabolite,
        MetaboliteWithCoefficient,
        Vector{MetaboliteWithCoefficient}
    },
)

Make a forward-only Reaction from substrates and products.

COBREXA.:↔Method
↔(
    substrates::Union{
        Nothing,
        Metabolite,
        MetaboliteWithCoefficient,
        Vector{MetaboliteWithCoefficient},
    },
    products::Union{
        Nothing,
        Metabolite,
        MetaboliteWithCoefficient,
        Vector{MetaboliteWithCoefficient}
    },
)

Make a bidirectional (reversible) Reaction from substrates and products.

COBREXA.add_gene!Method
add_gene!(model::StandardModel, genes::Gene)

Add gene to model based on gene id.

COBREXA.add_genes!Method
add_genes!(model::StandardModel, genes::Vector{Gene})

Add genes to model based on gene id.

COBREXA.add_metabolite!Method
add_metabolite!(model::StandardModel, met::Metabolite)

Add met to model based on metabolite id.

COBREXA.add_metabolites!Method
add_metabolites!(model::StandardModel, mets::Vector{Metabolite})

Add mets to model based on metabolite id.

COBREXA.add_reaction!Method
add_reaction!(model::StandardModel, rxn::Reaction)

Add rxn to model based on reaction id.

COBREXA.add_reactions!Method
add_reactions!(model::StandardModel, rxns::Vector{Reaction})

Add rxns to model based on reaction id.

COBREXA.remove_gene!Method
remove_gene!(
    model::StandardModel,
    id::Vector{String};
    knockout_reactions::Bool = false,
)

Remove gene with id from model. If knockout_reactions is true, then also constrain reactions that require the genes to function to carry zero flux.

Example

remove_gene!(model, "g1")
COBREXA.remove_genes!Method
remove_genes!(
    model::StandardModel,
    ids::Vector{String};
    knockout_reactions::Bool = false,
)

Remove all genes with ids from model. If knockout_reactions is true, then also constrain reactions that require the genes to function to carry zero flux.

Example

remove_genes!(model, ["g1", "g2"])
COBREXA.@add_reactions!Macro
@add_reactions!(model::Symbol, ex::Expr)

Shortcut to add multiple reactions and their lower and upper bounds

Call variants

@add_reactions! model begin
    reaction_name, reaction
end

@add_reactions! model begin
    reaction_name, reaction, lower_bound
end

@add_reactions! model begin
    reaction_name, reaction, lower_bound, upper_bound
end

Examples

@add_reactions! model begin
    "v1", nothing → A, 0, 500
    "v2", A ↔ B + C, -500
    "v3", B + C → nothing
end
COBREXA.add_model_with_exchangesMethod
add_model_with_exchanges(
    community::CoreModel,
    model::MetabolicModel,
    exchange_rxn_mets::Dict{String,String};
    model_name = "unknown_species",
    biomass_id = nothing,
)::CoreModel

Add model to community, which is a pre-existing community model with exchange reactions and metabolites in the dictionary exchange_rxn_mets. The model_name is appended to each reaction and metabolite, see join_with_exchanges. If biomass_id is specified then a biomass metabolite for model is also added to the resulting model. The column corresponding to the biomass_id reaction then produces this new biomass metabolite with unit coefficient. Note, exchange_rxn_ids and exchange_met_ids must already exist in the community model.

Example

community = add_model_with_exchanges(community,
    model,
    exchange_rxn_mets;
    model_name="species_2",
    biomass_id="BIOMASS_Ecoli_core_w_GAM")
COBREXA.add_objective!Method
add_objective!(
    community::CoreModel,
    objective_mets::Vector{String};
    objective_weights = Float64[],
    objective_column_index = 0,
)

Add an objective to the community model. Supply the string names of the objective metabolites in objective_mets. Optionally specify the weight to assign each metabolite in the objective function, if unassigned then equal weight is assumed. Also, optionally specify whether the objective already exists in the model by assigning objective_column_index. If unassigned then an objective column will be added, otherwise the column at objective_column_index will be updated.

Note, the weights are negated inside the function so that the objective metabolites are seen as reagents/substrates, not products in the reaction equation.

Example

add_objective!(model, ["met1", "met2"]) # adds a new column with weights = [1,1]
add_objective!(model, ["met1", "met2"]; objective_weights=[0.1, 0.9]) # adds a new column
add_objective!(model, ["met1", "met2"]; objective_weights=[0.1, 0.9], objective_column_index=10) # updates column 10
COBREXA.join_with_exchangesMethod
join_with_exchanges(models::Vector{M},
    exchange_rxn_mets::Dict{String, String};
    biomass_ids::Vector{String},
    model_names=String[]
)

Return a CoreModel representing the community model of models joined through their exchange reactions and metabolites in the dictionary exchange_rxn_mets, which maps exchange reactions to their associated metabolite. These exchange reactions and metabolites link model metabolites to environmental metabolites and reactions. Optionally specify model_names to append a specific name to each reaction and metabolite of an organism for easier reference (default is species_i for each model index i in models). Note, the bounds of the environmental variables are all set to zero. Thus, to run a simulation you need to constrain them appropriately. All the other bounds are inherited from the models used to construct the community model.

If biomass_ids is supplied, then a community model is returned that has an extra reaction added to the end of the stoichiometric matrix (last column) that can be assigned as the objective reaction. It also creates biomass "metabolites" that can be used in this objective reaction. Note, this reaction is unspecified, further action needs to be taken to specify it, e.g. assign weights to the last column of the stoichiometric matrix in the rows corresponding to the biomass metabolites.

To further clarify how this join works. Suppose you have 2 organisms with stoichiometric matrices S₁ and S₂ and you want to link them with exchange_rxn_mets = Dict(er₁ => em₁, er₂ => em₂, er₃ => em₃, ...). Then a new community stoichiometric matrix is constructed that looks like this:

            _      er₁  er₂  er₃  ...  b_
           |S₁                           |
           |   S₂                        |
        em₁|                             |
S   =   em₂|                             |
        em₃|                             |
        ...|                             |
        bm₁|                             |
        bm₂|_                           _|

The exchange reactions in each model get linked to environmental metabolites, emᵢ, and these get linked to environmental exchanges, erᵢ. These erᵢ behave like normal single organism exchange reactions. When biomass_ids are supplied, each model's biomass becomes a pseudo-metabolite (bmᵢ). These can be weighted in column b, called the community_biomass reaction in the community model, if desired. Refer to the tutorial if this is unclear.

Example

m1 = load_model(core_model_path)
m2 = load_model(CoreModel, core_model_path)

# need to list ALL the exchanges that will form part of the entire model
exchange_rxn_mets = Dict(k => first(keys(reaction_stoichiometry(m1, ex_rxn)))
    for filter(looks_like_exchange_reaction, reactions(m1)))

biomass_ids = ["BIOMASS_Ecoli_core_w_GAM", "BIOMASS_Ecoli_core_w_GAM"]

community = join_with_exchanges(
    [m1, m2],
    exchange_rxn_mets;
    biomass_ids = biomass_ids,
)

Model variant specifications

Analysis functions

COBREXA.envelope_latticeMethod
envelope_lattice(
    model::MetabolicModel,
    ridxs::Vector{Int};
    samples = 10,
    ranges = collect(zip(bounds(model)...))[ridxs],
    reaction_samples = fill(samples, length(ridxs)),
)

Create a lattice (list of "tick" vectors) for reactions at indexes ridxs in a model. Arguments samples, ranges, and reaction_samples may be optionally specified to customize the lattice creation process.

COBREXA.objective_envelopeMethod
objective_envelope(
    model::MetabolicModel,
    ridxs::Vector{Int},
    optimizer;
    lattice_args = (),
    lattice = envelope_lattice(model, ridxs; lattice_args...),
    kwargs...,
)

Compute an array of objective values for the model for rates of reactions specified ridxs fixed to a regular range of values between their respective lower and upper bounds.

This can be used to compute a "production envelope" of a metabolite; but generalizes to any specifiable objective and to multiple dimensions of the examined space. To retrieve a production envelope of any metabolite, set the objective coefficient vector of the model to a vector that contains a single 1 for the exchange reaction that "outputs" this metabolite, and run objective_envelope with the exchange reaction of the "parameter" metabolite specified in ridxs.

Returns a named tuple that contains lattice with reference values of the metabolites, and an N-dimensional array values with the computed objective values, where N is the number of specified reactions. Because of the increasing dimensionality, the computation gets very voluminous with increasing length of ridxs. The lattice for computing the optima can be supplied in the argument; by default it is created by envelope_lattice called on the model and reaction indexes. Additional arguments for the call to envelope_lattice can be optionally specified in lattice_args.

kwargs are internally forwarded to screen_optmodel_modifications.

Example

julia> m = load_model("test/downloaded/e_coli_core.xml");

julia> envelope = objective_envelope(m, ["R_EX_gln__L_e", "R_EX_fum_e"],
                                     Tulip.Optimizer;
                                     lattice_args=(samples=6,));

julia> envelope.lattice   # the reaction rates for which the optima were computed
2-element Vector{Vector{Float64}}:
 [0.0, 200.0, 400.0, 600.0, 800.0, 1000.0]
 [0.0, 200.0, 400.0, 600.0, 800.0, 1000.0]

julia> envelope.values   # the computed flux objective values for each reaction rate combination
6×6 Matrix{Float64}:
  0.873922   9.25815  17.4538  19.56   20.4121  20.4121
 13.0354    17.508    19.9369  21.894  22.6825  22.6825
 16.6666    18.6097   20.2847  21.894  22.6825  22.6825
 16.6666    18.6097   20.2847  21.894  22.6825  22.6825
 16.6666    18.6097   20.2847  21.894  22.6825  22.6825
 16.6666    18.6097   20.2847  21.894  22.6825  22.6825
COBREXA.flux_balance_analysisMethod
flux_balance_analysis(
    model::M,
    optimizer;
    modifications = [],
) where {M<:MetabolicModel}

Run flux balance analysis (FBA) on the model optionally specifying modifications to the problem. Basically, FBA solves this optimization problem:

max cᵀx
s.t. S x = b
     xₗ ≤ x ≤ xᵤ

See "Orth, J., Thiele, I. & Palsson, B. What is flux balance analysis?. Nat Biotechnol 28, 245-248 (2010). https://doi.org/10.1038/nbt.1614" for more information.

The optimizer must be set to a JuMP-compatible optimizer, such as GLPK.Optimizer or Tulip.Optimizer

Optionally, you may specify one or more modifications to be applied to the model before the analysis, such as change_optimizer_attribute, change_objective, and change_sense.

Returns an optimized JuMP model.

Example

model = load_model("e_coli_core.json")
solution = flux_balance_analysis(model, GLPK.optimizer)
value.(solution[:x])  # extract flux steady state from the optimizer

biomass_reaction_id = findfirst(model.reactions, "BIOMASS_Ecoli_core_w_GAM")

modified_solution = flux_balance_analysis(model, GLPK.optimizer;
    modifications=[change_objective(biomass_reaction_id)])
COBREXA.flux_balance_analysis_dictMethod
flux_balance_analysis_dict(model::MetabolicModel, args...)::Union{Dict{String, Float64},Nothing}

A variant of FBA that returns a dictionary assigning fluxes to reactions, if the solution is found. Arguments are passed to flux_balance_analysis.

COBREXA.flux_balance_analysis_vecMethod
flux_balance_analysis_vec(args...)::Union{Vector{Float64},Nothing}

A variant of FBA that returns a vector of fluxes in the same order as reactions of the model, if the solution is found.

Arguments are passed to flux_balance_analysis.

COBREXA._max_variability_fluxMethod
_max_variability_flux(opt_model, rid, ret)

Internal helper for maximizing reactions in optimization model.

COBREXA.flux_variability_analysisMethod
flux_variability_analysis(
    model::MetabolicModel,
    reactions::Vector{Int},
    optimizer;
    modifications = [],
    workers = [myid()],
    bounds = z -> (z,z),
    ret = objective_value,
)::Matrix{Float64}

Flux variability analysis solves a pair of optimization problems in model for each flux listed in reactions:

 min,max xᵢ
s.t. S x = b
    xₗ ≤ x ≤ xᵤ
     cᵀx ≥ bounds(Z₀)[1]
     cᵀx ≤ bounds(Z₀)[2]

where Z₀:= cᵀx₀ is the objective value of an optimal solution of the associated FBA problem (see flux_balance_analysis for a related analysis, also for explanation of the modifications argument).

The bounds is a user-supplied function that specifies the objective bounds for the variability optimizations, by default it restricts the flux objective value to the precise optimum reached in FBA. It can return -Inf and Inf in first and second pair to remove the limit. Use gamma_bounds and objective_bounds for simple bounds.

optimizer must be set to a JuMP-compatible optimizer. The computation of the individual optimization problems is transparently distributed to workers (see Distributed.workers()).

ret is a function used to extract results from optimized JuMP models of the individual reactions. By default, it calls and returns the value of JuMP.objective_value. More information can be extracted e.g. by setting it to a function that returns a more elaborate data structure; such as m -> (JuMP.objective_value(m), JuMP.value.(m[:x])).

Returns a matrix of extracted ret values for minima and maxima, of total size (length(reactions),2). The optimizer result status is checked with is_solved; nothing is returned if the optimization failed for any reason.

Example

model = load_model("e_coli_core.json")
flux_variability_analysis(model, [1, 2, 3, 42], GLPK.optimizer)
COBREXA.flux_variability_analysis_dictMethod
flux_variability_analysis_dict(
    model::MetabolicModel,
    optimizer;
    kwargs...
)

A variant of flux_variability_analysis that returns the individual maximized and minimized fluxes of all reactions as two dictionaries (of dictionaries). All keyword arguments except ret are passed through.

Example

mins, maxs = flux_variability_analysis_dict(
    model,
    Tulip.Optimizer;
    bounds = objective_bounds(0.99),
    modifications = [
        change_optimizer_attribute("IPM_IterationsLimit", 500),
        change_constraint("EX_glc__D_e"; lb = -10, ub = -10),
        change_constraint("EX_o2_e"; lb = 0, ub = 0),
    ],
)
COBREXA.parsimonious_flux_balance_analysisMethod
parsimonious_flux_balance_analysis(
    model::MetabolicModel,
    optimizer;
    modifications = [],
    qp_modifications = [],
    relax_bounds=[1.0, 0.999999, 0.99999, 0.9999, 0.999, 0.99],
)

Run parsimonious flux balance analysis (pFBA) on the model. In short, pFBA runs two consecutive optimization problems. The first is traditional FBA:

max cᵀx = μ
s.t. S x = b
     xₗ ≤ x ≤ xᵤ

And the second is a quadratic optimization problem:

min Σᵢ xᵢ²
s.t. S x = b
     xₗ ≤ x ≤ xᵤ
     μ = μ⁰

Where the optimal solution of the FBA problem, μ⁰, has been added as an additional constraint. See "Lewis, Nathan E, Hixson, Kim K, Conrad, Tom M, Lerman, Joshua A, Charusanti, Pep, Polpitiya, Ashoka D, Adkins, Joshua N, Schramm, Gunnar, Purvine, Samuel O, Lopez-Ferrer, Daniel, Weitz, Karl K, Eils, Roland, König, Rainer, Smith, Richard D, Palsson, Bernhard Ø, (2010) Omic data from evolved E. coli are consistent with computed optimal growth from genome-scale models. Molecular Systems Biology, 6. 390. doi: accession:10.1038/msb.2010.47" for more details.

pFBA gets the model optimum by standard FBA (using flux_balance_analysis with optimizer and modifications), then finds a minimal total flux through the model that still satisfies the (slightly relaxed) optimum. This is done using a quadratic problem optimizer. If the original optimizer does not support quadratic optimization, it can be changed using the callback in qp_modifications, which are applied after the FBA. See the documentation of fluxbalanceanalysis for usage examples of modifications.

Thhe optimum relaxation sequence can be specified in relax parameter, it defaults to multiplicative range of [1.0, 0.999999, ..., 0.99] of the original bound.

Returns an optimized model that contains the pFBA solution; or nothing if the optimization failed.

Example

model = load_model("e_coli_core.json")
optmodel = parsimonious_flux_balance_analysis(model, biomass, Gurobi.Optimizer)
value.(solution[:x])  # extract the flux from the optimizer
COBREXA._check_screen_mods_argsMethod
_check_screen_mods_args(mods, args, modsname)

Internal helper to check the presence and shape of modification and argument arrays in screen and pals.

COBREXA.screenMethod
function screen(
    model::MetabolicModel;
    variants::Maybe{Array{V,N}} = nothing,
    analysis,
    args::Maybe{Array{A,N}} = nothing,
    workers = [myid()],
)::Array where {V<:AbstractVector,A,N}

Take an array of model-modifying function vectors in variants, and execute the function analysis on all variants of the model specified by variants. The computation is distributed over worker IDs in workers. If args are supplied (as an array of the same size as the variants), they are forwarded as arguments to the corresponding analysis function calls.

The array of variants must contain vectors of single-parameter functions, these are applied to model in order. The functions must not modify the model, but rather return a modified copy. The copy should be made as shallow as possible, to increase memory efficiency of the process. Variant generators that modify the argument model in-place will cause unpredictable results. Refer to the definition of screen_variant for details.

The function analysis will receive a single argument (the modified model), together with arguments from args expanded by .... Supply an array of tuples or vectors to pass in multiple arguments to each function. If the argument values should be left intact (not expanded to multiple arguments), they must be wrapped in single-item tuple or vector.

The modification and analysis functions are transferred to workers as-is; all packages required to run them (e.g. the optimization solvers) must be loaded there. Typically, you want to use the macro @everywhere using MyFavoriteSolver from Distributed package for loading the solvers.

Return value

The results of running analysis are collected in to the resulting array, in a way that preserves the shape of the variants, similarly as with pmap.

The results of analysis function must be serializable, preferably made only from pure Julia structures, because they may be transferred over the network between the computation nodes. For that reason, functions that return whole JuMP models that contain pointers to allocated C structures (such as flux_balance_analysis used with GLPK or Gurobi otimizers) will generally not in this context.

Example

function reverse_reaction(i::Int)
    (model::CoreModel) -> begin
        mod = copy(model)
        mod.S[:,i] .*= -1   # this is unrealistic but sufficient for demonstration
        mod
    end
end

m = load_model(CoreModel, "e_coli_core.xml")

screen_variants(m,
           [
               [reverse_reaction(5)],
               [reverse_reaction(3), reverse_reaction(6)]
           ],
           mod -> mod.S[:,3])  # observe the changes in S

screen_variants(m,
    [
        [reverse_reaction(5)],
        [reverse_reaction(3), reverse_reaction(6)]
    ],
    mod -> flux_balance_analysis_vec(mod, GLPK.Optimizer))  # run analysis
COBREXA.screen_optmodel_modificationsMethod
function screen_optmodel_modifications(
    model::MetabolicModel,
    optimizer;
    common_modifications::V,
    modifications::Maybe{Array{V,N}} = nothing,
    args::Maybe{Array{T,N}} = nothing,
    analysis = screen_optimize_objective,
    workers = [myid()],
) where {V<:AbstractVector,T<:Tuple,N}

Screen multiple modifications of the same optimization model.

This function is potentially more efficient than screen because it avoids making variants of the model structure and remaking of the optimization model. On the other hand, modification functions need to keep the optimization model in a recoverable state (one that leaves the model usable for the next modification), which limits the possible spectrum of modifications applied.

Internally, model is distributed to workers and transformed into the optimization model using make_optimization_model. common_modifications are applied to the models at that point. Next, vectors of functions in modifications are consecutively applied, and the result of analysis function called on model are collected to an array of the same extent as modifications. Calls of analysis are optionally supplied with extra arguments from args expanded with ..., just like in screen.

Both the modification functions (in vectors) and the analysis function here have 2 base parameters (as opposed to 1 with screen): first is the model (carried through as-is), second is the prepared JuMP optimization model that may be modified and acted upon. As an example, you can use modification change_constraint and analysis screen_optimize_objective.

COBREXA.screen_variantFunction
screen_variant(model::MetabolicModel, variant::Vector, analysis, args = ())

Helper function for screen that applies all single-argument functions in variant to the model (in order from "first" to "last"), and executes analysis on the result.

Can be used to test model variants locally.

COBREXA.screen_variantsMethod
screen_variants(model, variants, analysis; workers=[myid()])

A shortcut for screen that only works with model variants.

Analysis modifications

COBREXA.add_crowding_constraintMethod
add_crowding_constraint(weights::Dict{Int64, Float64})

Adds a molecular crowding constraint to the optimization problem: ∑ wᵢ × vᵢ ≤ 1 where wᵢ is a weight and vᵢ is a flux index in the model's reactions specified in weights as vᵢ => wᵢ pairs.

See Beg, Qasim K., et al. "Intracellular crowding defines the mode and sequence of substrate uptake by Escherichia coli and constrains its metabolic activity." Proceedings of the National Academy of Sciences 104.31 (2007) for more details.

COBREXA.change_constraintMethod
change_constraint(id::String; lb=nothing, ub=nothing)

Change the lower and upper bounds (lb and ub respectively) of reaction id if supplied.

COBREXA.change_objectiveMethod
change_objective(new_objective::Union{String,Vector{String}}; weights=[], sense=MAX_SENSE)

Modification that changes the objective function used in a constraint based analysis function. new_objective can be a single reaction identifier, or an array of reactions identifiers.

Optionally, the objective can be weighted by a vector of weights, and a optimization sense can be set to either MAX_SENSE or MIN_SENSE.

COBREXA._do_knockoutMethod
_do_knockout(model::MetabolicModel, opt_model)

Internal helper for knockouts on generic MetabolicModels. This can be overloaded so that the knockouts may work differently (more efficiently) with other models.

COBREXA.knockoutMethod
knockout(gene_ids::Vector{String})

A modification that zeroes the bounds of all reactions that would be knocked out by the specified genes (effectively disables the reactions).

COBREXA.change_optimizerMethod
change_optimizer(optimizer)

Change the JuMP optimizer used to run the optimization.

This may be used to try different approaches for reaching the optimum, and in problems that may require different optimizers for different parts, such as the parsimonious_flux_balance_analysis.

COBREXA.change_optimizer_attributeMethod
change_optimizer_attribute(attribute_key, value)

Change a JuMP optimizer attribute. The attributes are optimizer-specific, refer to the JuMP documentation and the documentation of the specific optimizer for usable keys and values.

COBREXA.change_senseMethod
change_sense(objective_sense)

Change the objective sense of optimization. Possible arguments are MOI.MAX_SENSE and MOI.MIN_SENSE.

If you want to change the objective and sense at the same time, use change_objective instead to do both at once.

COBREXA.silenceFunction
silence

Modification that disable all output from the JuMP optimizer (shortcut for set_silent from JuMP).

Flux sampling

COBREXA._affine_hit_and_run_chainMethod
_affine_hit_and_run_chain(warmup, lbs, ubs, iters, chain)

Internal helper function for computing a single affine hit-and-run chain. The number of the chain is passed for possible future initialization of stable RNGs.

COBREXA.affine_hit_and_runMethod
function affine_hit_and_run(
    warmup_points::Matrix{Float64},
    lbs::Vector{Float64},
    ubs::Vector{Float64};
    sample_iters = 100 .* (1:5),
    workers = [myid()],
    chains = length(workers),
)

Run a hit-and-run style sampling that starts from warmup_points and uses their affine combinations for generating the run directions to sample the space delimited by lbs and ubs. The points that represent fluxes in warmup_points should be organized in columns, i.e. warmup_points[:,1] is the first warmup flux.

There are total chains of hit-and-run runs, each on a batch of size(warmup_points, 2) points. The runs are scheduled on workers, for good load balancing chains should be ideally much greater than length(workers).

Each run continues for maximum(sample_iters) iterations; the numbers in sample_iters represent the iterations at which the whole "current" batch of points is collected for output. For example, sample_iters=[1,4,5] causes the process run for 5 iterations, returning the sample batch that was produced by 1st, 4th and last (5th) iteration.

Returns a matrix of sampled fluxes (in columns), with all collected samples horizontally concatenated. The total number of samples (columns) will be size(warmup_points,2) * chains * length(sample_iters).

Example

using COBREXA
using Tulip

model = load_model(StandardModel, model_path)

warmup, lbs, ubs = warmup_from_variability(model, Tulip.Optimizer, 100)
samples = affine_hit_and_run(warmup, lbs, ubs, sample_iters = 1:3)
COBREXA.warmup_from_variabilityFunction
function warmup_from_variability(
    model::MetabolicModel,
    optimizer,
    min_reactions::Vector{Int}=1:n_reactions(model),
    max_reactions::Vector{Int}=1:n_reactions(model);
    modifications = [],
    workers::Vector{Int} = [myid()],
)::Tuple{Matrix{Float64}, Vector{Float64}, Vector{Float64}}

Generate FVA-like warmup points for samplers, by minimizing and maximizing the specified reactions. The result is returned as a matrix, each point occupies as single column in the result.

COBREXA.warmup_from_variabilityMethod
warmup_from_variability(
    model::MetabolicModel,
    optimizer,
    n_points::Int;
    kwargs...
)

Generates FVA-like warmup points for samplers, by selecting random points by minimizing and maximizing reactions. Can not return more than 2 times the number of reactions in the model.

Miscellaneous utilities

COBREXA.ambiguously_identified_itemsMethod
ambiguously_identified_items(
    index::Dict{String,Dict{String,[String]}},
)::Vector{String}

Find items (genes, metabolites, ...) from the annotation index that are identified non-uniquely by at least one of their annotations.

This often indicates that the items are duplicate or miscategorized.

COBREXA.annotation_indexMethod
annotation_index(
    xs::AbstractDict{String};
    annotations = _annotations,
)::Dict{String,Dict{String,[String]}}

Extract annotations from a dictionary of items xs and build an index that maps annotation "kinds" (e.g. "PubChem") to the mapping from the annotations (e.g. "COMPOUND_12345") to item IDs that carry the annotations.

Function annotations is used to access the Annotations object in the dictionary values.

This is extremely useful for finding items by annotation data.

COBREXA.check_duplicate_reactionMethod
check_duplicate_reaction(rxn::Reaction, rxns::Dict{String, Reaction}; only_metabolites=true)

Check if rxn already exists in rxns but has another id. If only_metabolites is true then only the metabolite ids are checked. Otherwise, compares metabolite ids and the absolute value of their stoichiometric coefficients to those of rxn. If rxn has the same reaction equation as another reaction in rxns, the return the id. Otherwise return nothing.

See also: reaction_mass_balanced

COBREXA.is_boundaryMethod
is_boundary(reaction::Reaction)

Return true if reaction is a boundary reaction, otherwise return false. Checks if on boundary by inspecting the number of metabolites in reaction equation. Boundary reactions have only one metabolite, e.g. an exchange reaction, or a sink/demand reaction.

COBREXA.reaction_atom_balanceMethod
reaction_atom_balance(model::StandardModel, rxn)

Returns a dictionary mapping the stoichiometry of atoms through a single reaction. Uses the metabolite information in model to determine the mass balance. Accepts a reaction dictionary, a reaction string id or a Reaction as an argument for rxn.

See also: reaction_mass_balanced

COBREXA.stoichiometry_stringMethod
stoichiometry_string(rxn_dict::Dict{String, Float64}; format_id = x -> x)

Return the reaction equation as a string. The metabolite strings can be manipulated by setting format_id.

Example

julia> req = Dict("coa_c" => -1, "for_c" => 1, "accoa_c" => 1, "pyr_c" => -1)
julia> stoichiometry_string(req)
"coa_c + pyr_c = for_c + accoa_c"

julia> stoichiometry_string(req; format_id = x -> x[1:end-2])
"coa + pyr = for + accoa"
COBREXA.serialize_modelMethod
serialize_model(model::Serialized, filename::String)::Serialized

Specialization of serialize_model that prevents nested serialization of already-serialized models.

COBREXA.serialize_modelMethod
serialize_model(model::MM, filename::String)::Serialized{MM} where {MM<:MetabolicModel}

Serialize the model to file filename, returning a Serialized model that is able to load itself back automatically upon precaching by precache!.

COBREXA.gamma_boundsMethod
gamma_bounds(gamma)

A bounds-generating function for flux_variability_analysis that limits the objective value to be at least gamma*Z₀, as usual in COBRA packages. Use as the bounds argument:

flux_variability_analysis(model, some_optimizer; bounds = gamma_bounds(0.9))
COBREXA._parse_formulaMethod
_parse_formula(f::String)::MetaboliteFormula

Parse a formula in format C2H6O into a MetaboliteFormula, which is basically a dictionary of atom counts in the molecule.

COBREXA.atom_fluxesMethod
atom_fluxes(model::MetabolicModel, reaction_fluxes::Dict{String, Float64})

Return a dictionary mapping the flux of atoms across a flux solution given by reaction_fluxes using the reactions in model to determine the appropriate stoichiometry.

Note, this function ignores metabolites with no formula assigned to them, no error message will be displayed.

Note, if a model is mass balanced there should be not net flux of any atom. By removing reactions from the flux_solution you can investigate how that impacts the mass balances.

Example

# Find flux of Carbon through all metabolic reactions except the biomass reaction
delete!(fluxes, "BIOMASS_Ecoli_core_w_GAM")
atom_fluxes(model, fluxes)["C"]
COBREXA.metabolite_fluxesMethod
metabolite_fluxes(model::MetabolicModel, flux_dict::Dict{String, Float64})

Return two dictionaries of metabolite ids mapped to reactions that consume or produce them, given the flux distribution supplied in flux_dict.

COBREXA._parse_grrMethod
_parse_grr(gpa::SBML.GeneProductAssociation)::GeneAssociation

Parse SBML.GeneProductAssociation structure to the simpler GeneAssociation. The input must be (implicitly) in a positive DNF.

COBREXA._parse_grrMethod
_parse_grr(s::String)::GeneAssociation

Parse a DNF gene association rule in format (YIL010W and YLR043C) or (YIL010W and YGR209C) to GeneAssociation. Also acceptsOR,|,||,AND,&, and&&`.

Example

julia> _parse_grr("(YIL010W and YLR043C) or (YIL010W and YGR209C)")
2-element Array{Array{String,1},1}:
 ["YIL010W", "YLR043C"]
 ["YIL010W", "YGR209C"]
COBREXA._unparse_grrMethod
_unparse_grr(
    ::Type{SBML.GeneProductAssociation},
    x::GeneAssociation,
)::SBML.GeneAssociation

Convert a GeneAssociation to the corresponding SBML.jl structure.

COBREXA._unparse_grrMethod
unparse_grr(grr::Vector{Vector{Gene}}

Converts a nested string gene reaction array back into a gene reaction rule string.

Example

julia> _unparse_grr(String, [["YIL010W", "YLR043C"], ["YIL010W", "YGR209C"]])
"(YIL010W and YLR043C) or (YIL010W and YGR209C)"
COBREXA._guesskeyMethod
_guesskey(ks, possibilities)

Unfortunately, many model types that contain dictionares do not have standardized field names, so we need to try a few possibilities and guess the best one. The keys used to look for valid field names should be ideally specified as constants in src/base/constants.jl.

COBREXA.looks_like_biomass_reactionMethod
looks_like_biomass_reaction(rxn_id::String;
    exclude_exchanges = false,
    exchange_prefixes = _constants.exchange_prefixes,
    biomass_strings = _constants.biomass_strings,
)::Bool

A predicate that matches reaction identifiers that look like biomass reactions. Biomass reactions are identified by looking for occurences of biomass_strings in the reaction id. If exclude_exchanges is set, the strings that look like exchanges (from looks_like_exchange_reaction) will not match.

Example

filter(looks_like_biomass_reaction, reactions(model)) # returns strings
findall(looks_like_biomass_reaction, reactions(model)) # returns indices
COBREXA.looks_like_exchange_metaboliteMethod
looks_like_exchange_metabolite(rxn_id::String;
    exchange_suffixes = _constants.exchange_suffixes,
    )::Bool

A predicate that matches metabolite identifiers that look like involved in exchange reactions. Exchange metabolites are identified by exchange_suffixes at the end of the metabolite id.

Example

filter(looks_like_exchange_metabolite, metabolites(model)) # returns strings
findall(looks_like_exchange_metabolite, metabolites(model)) # returns indices
COBREXA.looks_like_exchange_reactionMethod
looks_like_exchange_reaction(rxn_id::String;
    exclude_biomass = false,
    biomass_strings = _constants.biomass_strings,
    exchange_prefixes = _constants.exchange_prefixes,
)

A predicate that matches reaction identifiers that look like exchange or biomass reactions, given the usual naming schemes in common model repositories. Exchange reactions are identified based on matching prefixes in the set exchange_prefixes and biomass reactions are identified by looking for occurences of biomass_strings in the reaction id.

Also see find_exchange_reactions.

Example

findall(looks_like_exchange_reaction, reactions(model)) # returns indices
filter(looks_like_exchange_reaction, reactions(model)) # returns Strings

# to use the optional arguments you need to expand the function's arguments
# using an anonymous function
findall(x -> looks_like_exchange_reaction(x; exclude_biomass=true), reactions(model)) # returns indices
filter(x -> looks_like_exchange_reaction(x; exclude_biomass=true), reactions(model)) # returns Strings
COBREXA.looks_like_internal_reactionMethod
looks_like_internal_reaction(
    rxn_id::String;
    exchange_prefixes = _constants.exchange_prefixes,
    biomass_strings = _constants.biomass_strings,
)::Bool

A predicate that matches reaction identifiers that look like internal reactions, i.e. reactions that are neither exchange nor biomass reactions. Exchange reactions are identified based on matching prefixes in the set exchange_prefixes and biomass reactions are identified by looking for occurences of biomass_strings in the reaction id.

Also see find_internal_reactions.

Example

findall(looks_like_internal_reaction, reactions(model)) # returns indices
filter(looks_like_internal_reaction, reactions(model)) # returns Strings

# to use the optional arguments you need to expand the function's arguments
# using an anonymous function
findall(x -> looks_like_internal_reaction(x; exchange_prefixes=["EX_", "R_EX_"]), reactions(model)) # returns indices
filter(x -> looks_like_internal_reaction(x; exchange_prefixes=["EX_", "R_EX_"]), reactions(model)) # returns Strings
COBREXA.change_bound!Method
change_bound!(
    model::CoreModel,
    rxn_idx::Int;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model in-place.

Example

new_model = change_bound!(model, 123, lower=-21.15, upper=42.3)
COBREXA.change_bound!Method
change_bound!(
    model::CoreModel,
    rxn_id::String;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model in-place.

Example

new_model = change_bound!(model, "ReactionB", lower=-21.15, upper=42.3)
COBREXA.change_bound!Method
change_bound!(
    model::CoreModelCoupled,
    rxn_idx::Int;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model in-place.

Example

new_model = change_bound!(model, 123, lower=-21.15, upper=42.3)
COBREXA.change_bound!Method
change_bound!(
    model::CoreModelCoupled,
    rxn_id::String;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model in-place.

Example

new_model = change_bound!(model, "ReactionB", lower=-21.15, upper=42.3)
COBREXA.change_bound!Method
change_bound!(
    model::StandardModel,
    rxn_id::String;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model in-place.

Example

new_model = change_bound!(model, "ReactionB", lower=-21.15, upper=42.3)
COBREXA.change_boundMethod
change_bound(
    model::CoreModel,
    rxn_idx::Int;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model and return the modified model.

Example

change_bound(model, 123, lower=-21.15, upper=42.3)
COBREXA.change_boundMethod
change_bound(
    model::CoreModel,
    rxn_id::String;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model and return the modified model.

Example

change_bound(model, "ReactionB", lower=-21.15, upper=42.3)
COBREXA.change_boundMethod
change_bound(
    model::CoreModelCoupled,
    rxn_idx::Int;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model and return the modified model.

Example

change_bound(model, 123, lower=-21.15, upper=42.3)
COBREXA.change_boundMethod
change_bound(
    model::CoreModelCoupled,
    rxn_id::String;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model and return the modified model.

Example

change_bound(model, "ReactionB", lower=-21.15, upper=42.3)
COBREXA.change_boundMethod
change_bound(
    model::StandardModel,
    rxn_id::String;
    lower = nothing,
    upper = nothing,
)

Change the specified reaction flux bound in the model and return the modified model.

Example

change_bound(model, "ReactionB", lower=-21.15, upper=42.3)
COBREXA.change_bounds!Method
change_bounds!(
    model::CoreModel,
    rxn_idxs::AbstractVector{Int64};
    lower = (nothing for _ = rxn_idxs),
    upper = (nothing for _ = rxn_idxs),
)

Change the specified reaction flux bounds in the model in-place.

Example

new_model = change_bounds!(model, [123, 234], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_bounds!Method
change_bounds!(
    model::CoreModel,
    rxn_ids::AbstractVector{String};
    lower = (nothing for _ = rxn_ids),
    upper = (nothing for _ = rxn_ids),
)

Change the specified reaction flux bounds in the model in-place.

Example

new_model = change_bounds!(model, ["ReactionA", "ReactionC"], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_bounds!Method
change_bounds!(
    model::CoreModelCoupled,
    rxn_idxs::AbstractVector{Int64};
    lower = (nothing for _ = rxn_idxs),
    upper = (nothing for _ = rxn_idxs),
)

Change the specified reaction flux bounds in the model in-place.

Example

new_model = change_bounds!(model, [123, 234], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_bounds!Method
change_bounds!(
    model::CoreModelCoupled,
    rxn_ids::AbstractVector{String};
    lower = (nothing for _ = rxn_ids),
    upper = (nothing for _ = rxn_ids),
)

Change the specified reaction flux bounds in the model in-place.

Example

new_model = change_bounds!(model, ["ReactionA", "ReactionC"], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_bounds!Method
change_bounds!(
    model::StandardModel,
    rxn_ids::AbstractVector{String};
    lower = (nothing for _ = rxn_ids),
    upper = (nothing for _ = rxn_ids),
)

Change the specified reaction flux bounds in the model in-place.

Example

new_model = change_bounds!(model, ["ReactionA", "ReactionC"], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_boundsMethod
change_bounds(
    model::CoreModel,
    rxn_idxs::AbstractVector{Int64};
    lower = (nothing for _ = rxn_idxs),
    upper = (nothing for _ = rxn_idxs),
)

Change the specified reaction flux bounds in the model and return the modified model.

Example

change_bounds(model, [123, 234], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_boundsMethod
change_bounds(
    model::CoreModel,
    rxn_ids::AbstractVector{String};
    lower = (nothing for _ = rxn_ids),
    upper = (nothing for _ = rxn_ids),
)

Change the specified reaction flux bounds in the model and return the modified model.

Example

change_bounds(model, ["ReactionA", "ReactionC"], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_boundsMethod
change_bounds(
    model::CoreModelCoupled,
    rxn_idxs::AbstractVector{Int64};
    lower = (nothing for _ = rxn_idxs),
    upper = (nothing for _ = rxn_idxs),
)

Change the specified reaction flux bounds in the model and return the modified model.

Example

change_bounds(model, [123, 234], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_boundsMethod
change_bounds(
    model::CoreModelCoupled,
    rxn_ids::AbstractVector{String};
    lower = (nothing for _ = rxn_ids),
    upper = (nothing for _ = rxn_ids),
)

Change the specified reaction flux bounds in the model and return the modified model.

Example

change_bounds(model, ["ReactionA", "ReactionC"], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.change_boundsMethod
change_bounds(
    model::StandardModel,
    rxn_ids::AbstractVector{String};
    lower = (nothing for _ = rxn_ids),
    upper = (nothing for _ = rxn_ids),
)

Change the specified reaction flux bounds in the model and return the modified model.

Example

change_bounds(model, ["ReactionA", "ReactionC"], lower=[-2.1, -50.05], upper=[4.2, 100.1])
COBREXA.remove_metabolite!Method
remove_metabolite!(model::CoreModel, metabolite_idx::Int)

Remove metabolite from the model of type CoreModel in-place.

COBREXA.remove_metabolite!Method
remove_metabolite!(model::CoreModel, metabolite_id::String)

Remove metabolite from the model of type CoreModel in-place.

COBREXA.remove_metabolite!Method
remove_metabolite!(model::CoreModelCoupled, metabolite_idx::Int)

Remove metabolite from the model of type CoreModelCoupled in-place.

COBREXA.remove_metabolite!Method
remove_metabolite!(model::CoreModelCoupled, metabolite_id::String)

Remove metabolite from the model of type CoreModelCoupled in-place.

COBREXA.remove_metabolite!Method
remove_metabolite!(model::StandardModel, metabolite_id::String)

Remove metabolite from the model of type StandardModel in-place.

COBREXA.remove_metaboliteMethod
remove_metabolite(model::CoreModel, metabolite_idx::Int)

Remove metabolite from the model of type CoreModel and return the modified model.

COBREXA.remove_metaboliteMethod
remove_metabolite(model::CoreModel, metabolite_id::String)

Remove metabolite from the model of type CoreModel and return the modified model.

COBREXA.remove_metaboliteMethod
remove_metabolite(model::CoreModelCoupled, metabolite_idx::Int)

Remove metabolite from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_metaboliteMethod
remove_metabolite(model::CoreModelCoupled, metabolite_id::String)

Remove metabolite from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_metaboliteMethod
remove_metabolite(model::StandardModel, metabolite_id::String)

Remove metabolite from the model of type StandardModel and return the modified model.

COBREXA.remove_metabolites!Method
remove_metabolites!(model::CoreModel, metabolite_idxs::AbstractVector{Int64})

Remove metabolites from the model of type CoreModel in-place.

COBREXA.remove_metabolites!Method
remove_metabolites!(model::CoreModel, metabolite_ids::AbstractVector{String})

Remove metabolites from the model of type CoreModel in-place.

COBREXA.remove_metabolites!Method
remove_metabolites!(model::CoreModelCoupled, metabolite_idxs::AbstractVector{Int64})

Remove metabolites from the model of type CoreModelCoupled in-place.

COBREXA.remove_metabolites!Method
remove_metabolites!(model::CoreModelCoupled, metabolite_ids::AbstractVector{String})

Remove metabolites from the model of type CoreModelCoupled in-place.

COBREXA.remove_metabolites!Method
remove_metabolites!(model::StandardModel, metabolite_ids::AbstractVector{String})

Remove metabolites from the model of type StandardModel in-place.

COBREXA.remove_metabolitesMethod
remove_metabolites(model::CoreModel, metabolite_idxs::AbstractVector{Int64})

Remove metabolites from the model of type CoreModel and return the modified model.

COBREXA.remove_metabolitesMethod
remove_metabolites(model::CoreModel, metabolite_ids::AbstractVector{String})

Remove metabolites from the model of type CoreModel and return the modified model.

COBREXA.remove_metabolitesMethod
remove_metabolites(model::CoreModelCoupled, metabolite_idxs::AbstractVector{Int64})

Remove metabolites from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_metabolitesMethod
remove_metabolites(model::CoreModelCoupled, metabolite_ids::AbstractVector{String})

Remove metabolites from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_metabolitesMethod
remove_metabolites(model::StandardModel, metabolite_ids::AbstractVector{String})

Remove metabolites from the model of type StandardModel and return the modified model.

COBREXA.remove_reaction!Method
remove_reaction!(model::CoreModel, reaction_idx::Int)

Remove reaction from the model of type CoreModel in-place.

COBREXA.remove_reaction!Method
remove_reaction!(model::CoreModel, reaction_id::String)

Remove reaction from the model of type CoreModel in-place.

COBREXA.remove_reaction!Method
remove_reaction!(model::CoreModelCoupled, reaction_idx::Int)

Remove reaction from the model of type CoreModelCoupled in-place.

COBREXA.remove_reaction!Method
remove_reaction!(model::CoreModelCoupled, reaction_id::String)

Remove reaction from the model of type CoreModelCoupled in-place.

COBREXA.remove_reaction!Method
remove_reaction!(model::StandardModel, reaction_id::String)

Remove reaction from the model of type StandardModel in-place.

COBREXA.remove_reactionMethod
remove_reaction(model::CoreModel, reaction_idx::Int)

Remove reaction from the model of type CoreModel and return the modified model.

COBREXA.remove_reactionMethod
remove_reaction(model::CoreModel, reaction_id::String)

Remove reaction from the model of type CoreModel and return the modified model.

COBREXA.remove_reactionMethod
remove_reaction(model::CoreModelCoupled, reaction_idx::Int)

Remove reaction from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_reactionMethod
remove_reaction(model::CoreModelCoupled, reaction_id::String)

Remove reaction from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_reactionMethod
remove_reaction(model::StandardModel, reaction_id::String)

Remove reaction from the model of type StandardModel and return the modified model.

COBREXA.remove_reactions!Method
remove_reactions!(model::CoreModel, reaction_idxs::AbstractVector{Int64})

Remove reactions from the model of type CoreModel in-place.

COBREXA.remove_reactions!Method
remove_reactions!(model::CoreModel, reaction_ids::AbstractVector{String})

Remove reactions from the model of type CoreModel in-place.

COBREXA.remove_reactions!Method
remove_reactions!(model::CoreModelCoupled, reaction_idxs::AbstractVector{Int64})

Remove reactions from the model of type CoreModelCoupled in-place.

COBREXA.remove_reactions!Method
remove_reactions!(model::CoreModelCoupled, reaction_ids::AbstractVector{String})

Remove reactions from the model of type CoreModelCoupled in-place.

COBREXA.remove_reactions!Method
remove_reactions!(model::StandardModel, reaction_ids::AbstractVector{String})

Remove reactions from the model of type StandardModel in-place.

COBREXA.remove_reactionsMethod
remove_reactions(model::CoreModel, reaction_idxs::AbstractVector{Int64})

Remove reactions from the model of type CoreModel and return the modified model.

COBREXA.remove_reactionsMethod
remove_reactions(model::CoreModel, reaction_ids::AbstractVector{String})

Remove reactions from the model of type CoreModel and return the modified model.

COBREXA.remove_reactionsMethod
remove_reactions(model::CoreModelCoupled, reaction_idxs::AbstractVector{Int64})

Remove reactions from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_reactionsMethod
remove_reactions(model::CoreModelCoupled, reaction_ids::AbstractVector{String})

Remove reactions from the model of type CoreModelCoupled and return the modified model.

COBREXA.remove_reactionsMethod
remove_reactions(model::StandardModel, reaction_ids::AbstractVector{String})

Remove reactions from the model of type StandardModel and return the modified model.

Logging and debugging helpers

COBREXA.log_ioFunction
log_io(enable::Bool=true)

Enable (default) or disable (by passing false) output of messages and warnings from model input/output.

COBREXA.log_modelsFunction
log_models(enable::Bool=true)

Enable (default) or disable (by passing false) output of model-related messages.

COBREXA.log_perfFunction
log_perf(enable::Bool=true)

Enable (default) or disable (by passing false) output of performance-related tracing information.

COBREXA.@_make_logging_tagMacro
macro _make_logging_group(sym::Symbol, doc::String)

This creates a group of functions that allow masking out topic-related logging actions. A call that goes as follows:

@_make_logging_tag XYZ

creates the following tools:

  • global variable _XYZ_log_enabled defaulted to false
  • function log_XYZ that can be called to turn the logging on/off
  • a masking macro @_XYZ_log that can be prepended to commands that should only happen if the logging of tag XYZ is enabled.

The masking macro is then used as follows:

@_XYZ_log @info "This is the extra verbose information you wanted!" a b c

The user can direct logging with these:

log_XYZ()
log_XYZ(false)

doc should be a name of the stuff that is being printed if the corresponding log_XYZ() is enabled – it is used to create a friendly documentation for the logging switch. In this case it could say "X, Y and Z-related messages".