Utilities

Helper functions

COBREXA.ambiguously_identified_itemsMethod
ambiguously_identified_items(index::Dict{String, Dict{String, Set{String}}}) -> Set{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.

source
COBREXA.annotation_indexMethod
annotation_index(xs::AbstractDict{String}; annotations) -> Dict{String, Dict{String, Set{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.

source
COBREXA.check_duplicate_reactionMethod
check_duplicate_reaction(crxn::Reaction, rxns::OrderedCollections.OrderedDict{String, Reaction}; only_metabolites) -> Union{Nothing, String}

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

source
COBREXA.is_boundaryMethod
is_boundary(rxn_dict::Dict{String, Float64}) -> Bool

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

source
COBREXA.reaction_atom_balanceMethod
reaction_atom_balance(model::MetabolicModel, reaction_dict::Dict{String, Float64}) -> Dict{String, Float64}

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

source
COBREXA.stoichiometry_stringMethod
stoichiometry_string(req; format_id) -> String

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"
source
COBREXA.serialize_modelMethod
serialize_model(model::MetabolicModel, filename::String) -> Serialized

Serialize the model to file filename, returning a Serialized model that can be loaded back transparently by precache!. The result does not contain the actual model data that are deferred to the disk; it may thus be used to save memory, or send the model efficiently to remote workers within distributed shared-storage environments.

The benefit of using this over "raw" Serialization.serialize is that the resulting Serialized model will reload itself automatically with precache! at first usage, which needs to be done manually when using the Serialization package directly.

source
COBREXA.gamma_boundsMethod
gamma_bounds(gamma) -> COBREXA.var"#659#660"

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))
source
COBREXA.expression_probabilistic_boundsMethod
expression_probabilistic_bounds(relative_expression::Dict{String, Float64}, grr::Union{Nothing, Vector{Vector{String}}}) -> Float64

Create E-Flux-like bounds for a reaction that requires gene products given by grr, with expression "choke" data given by relative probabilities (between 0 and 1) in relative_expression.

source
COBREXA.gene_product_dictMethod
gene_product_dict(model::GeckoModel, opt_model) -> Any

Return a dictionary mapping protein molar concentrations to their ids. The argument opt_model is a solved optimization problem, typically returned by flux_balance_analysis. See flux_dict for the corresponding function that returns a dictionary of solved fluxes.

source
COBREXA.atom_fluxesMethod
atom_fluxes(model::MetabolicModel, reaction_fluxes::Dict{String, Float64}) -> 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"]
source
COBREXA.metabolite_fluxesMethod
metabolite_fluxes(model::MetabolicModel, flux_dict::Dict{String, Float64}) -> Tuple{Dict{String, Dict{String, Float64}}, Dict{String, 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.

source
COBREXA._gecko_gene_product_couplingMethod
_gecko_gene_product_coupling(model::GeckoModel) -> SparseArrays.SparseMatrixCSC{Float64, Int64}

Compute the part of the coupling for GeckoModel that limits the amount of each kind of protein available.

source
COBREXA._gecko_mass_group_couplingMethod
_gecko_mass_group_coupling(model::GeckoModel) -> SparseArrays.SparseMatrixCSC{Float64, Int64}

Compute the part of the coupling for GeckoModel that limits the total mass of each group of gene products.

source
COBREXA._gecko_reaction_column_reactionsMethod
_gecko_reaction_column_reactions(model::GeckoModel) -> SparseArrays.SparseMatrixCSC{Int64, Int64}

Retrieve a utility mapping between reactions and split reactions; rows correspond to "original" reactions, columns correspond to "split" reactions.

source
COBREXA._gecko_reaction_couplingMethod
_gecko_reaction_coupling(model::GeckoModel) -> SparseArrays.SparseMatrixCSC{Int64, Int64}

Compute the part of the coupling for GeckoModel that limits the "arm" reactions (which group the individual split unidirectional reactions).

source
COBREXA._parse_grrMethod
_parse_grr(gpa::SBML.GeneProductAssociation) -> Vector{Vector{String}}

Parse SBML.GeneProductAssociation structure and convert it to a strictly positive DNF GeneAssociation. Negation (SBML.GPANot) is not supported.

source
COBREXA._parse_grrMethod
_parse_grr(s::String) -> Union{Nothing, Vector{Vector{String}}}

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"]
source
COBREXA._parse_grr_to_sbmlMethod
_parse_grr_to_sbml(str::String) -> Union{Nothing, SBML.GeneProductAssociation}

Internal helper for parsing the string GRRs into SBML data structures. More general than _parse_grr.

source
COBREXA._sortuniqueMethod
_sortunique(x) -> Any

A helper for producing predictable unique sequences. Might be faster if compacting would be done directly in sort().

source
COBREXA._unparse_grrMethod
_unparse_grr(_::Type{SBML.GeneProductAssociation}, x::Vector{Vector{String}}) -> SBML.GPAOr

Convert a GeneAssociation to the corresponding SBML.jl structure.

source
COBREXA._unparse_grrMethod
_unparse_grr(::Type{String}, grr::Vector{Vector{String}}; and, or) -> String

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

Example

julia> _unparse_grr(String, [["YIL010W", "YLR043C"], ["YIL010W", "YGR209C"]])
"(YIL010W && YLR043C) || (YIL010W && YGR209C)"
source
COBREXA._guesskeyMethod
_guesskey(avail, possibilities) -> Any

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.

source
COBREXA.getsMethod
gets(collection, fail, keys) -> Any

Return fail if key in keys is not in collection, otherwise return collection[key]. Useful if may different keys need to be tried due to non-standardized model formats.

source
COBREXA.is_atp_maintenance_reactionMethod
is_atp_maintenance_reaction(
    model::MetabolicModel,
    reaction_id::String;
    annotation_keys = ["sbo", "SBO"],
)

Check if a reaction is annotated as an atpmaintenance reaction. Uses `Identifiers.ATPMAINTENANCEREACTIONSinternally, which includes SBO identifiers. In the reaction annotations, use the keys inannotationkeys` to look for entries. Returns false if no hits or if no keys are found.

source
COBREXA.is_biomass_reactionMethod
is_biomass_reaction(
    model::MetabolicModel,
    reaction_id::String;
    annotation_keys = ["sbo", "SBO"],
)

Check if a reaction is annotated as a biomass reaction. Uses Identifiers.BIOMASS_REACTIONS internally, which includes SBO identifiers. In the reaction annotations, use the keys in annotation_keys to look for entries. Returns false if no hits or if no keys are found.

source
COBREXA.is_exchange_reactionMethod
is_exchange_reaction(
    model::MetabolicModel,
    reaction_id::String;
    annotation_keys = ["sbo", "SBO"],
)

Check if a reaction is annotated as an exchange reaction. Uses Identifiers.EXCHANGE_REACTIONS internally, which includes SBO identifiers. In the reaction annotations, use the keys in annotation_keys to look for entries. Returns false if no hits or if no keys are found.

source
COBREXA.is_metabolic_reactionMethod
is_metabolic_reaction(
    model::MetabolicModel,
    reaction_id::String;
    annotation_keys = ["sbo", "SBO"],
)

Check if a reaction is annotated as a metabolic reaction. Uses Identifiers.METABOLIC_REACTIONS internally, which includes SBO identifiers. In the reaction annotations, use the keys in annotation_keys to look for entries. Returns false if no hits or if no keys are found.

source
COBREXA.is_pseudo_reactionMethod
is_pseudo_reaction(
    model::MetabolicModel,
    reaction_id::String;
    annotation_keys = ["sbo", "SBO"],
)

Check if a reaction is annotated as a pseudo reaction. Uses Identifiers.PSEUDOREACTIONS internally, which includes SBO identifiers. In the reaction annotations, use the keys in annotation_keys to look for entries. Returns false if no hits or if no keys are found.

source
COBREXA.is_spontaneous_reactionMethod
is_spontaneous_reaction(
    model::MetabolicModel,
    reaction_id::String;
    annotation_keys = ["sbo", "SBO"],
)

Check if a reaction is annotated as a spontaneous reaction. Uses Identifiers.SPONTANEOUS_REACTIONS internally, which includes SBO identifiers. In the reaction annotations, use the keys in annotation_keys to look for entries. Returns false if no hits or if no keys are found.

source
COBREXA.is_transport_reactionMethod
is_transport_reaction(
    model::MetabolicModel,
    reaction_id::String;
    annotation_keys = ["sbo", "SBO"],
)

Check if a reaction is annotated as a transport reaction. Uses Identifiers.TRANSPORT_REACTIONS internally, which includes SBO identifiers. In the reaction annotations, use the keys in annotation_keys to look for entries. Returns false if no hits or if no keys are found.

source
COBREXA.looks_like_biomass_reactionMethod
looks_like_biomass_reaction(rxn_id::String; exclude_exchanges, exchange_prefixes, 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.

Note

While looks_like_biomass_reaction is useful for heuristically finding a reaction, it is preferable to use standardized terms for finding reactions (e.g. SBO terms). See is_biomass_reaction for a more systematic alternative.

Example

filter(looks_like_biomass_reaction, reactions(model)) # returns strings
findall(looks_like_biomass_reaction, reactions(model)) # returns indices
source
COBREXA.looks_like_exchange_reactionMethod
looks_like_exchange_reaction(rxn_id::String; exclude_biomass, biomass_strings, exchange_prefixes) -> Bool

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.

Note

While looks_like_exchange_reaction is useful for heuristically finding a reaction, it is preferable to use standardized terms for finding reactions (e.g. SBO terms). See is_exchange_reaction for a more systematic alternative.

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
source
COBREXA.looks_like_extracellular_metaboliteMethod
looks_like_extracellular_metabolite(met_id::String; extracellular_suffixes) -> Bool

A predicate that matches metabolite identifiers that look like they are extracellular metabolites. Extracellular metabolites are identified by extracellular_suffixes at the end of the metabolite id.

Example

filter(looks_like_extracellular_metabolite, metabolites(model)) # returns strings
findall(looks_like_extracellular_metabolite, metabolites(model)) # returns indices
source
COBREXA._smoment_column_reactionsMethod
_smoment_column_reactions(model::SMomentModel) -> SparseArrays.SparseMatrixCSC{Int64, Int64}

Retrieve a utility mapping between reactions and split reactions; rows correspond to "original" reactions, columns correspond to "split" reactions.

source
COBREXA.smoment_isozyme_speedMethod
smoment_isozyme_speed(gene_product_molar_mass::Function) -> COBREXA.var"#775#776"

A piping- and argmax-friendly overload of smoment_isozyme_speed.

Example

gene_mass_function = gid -> 1.234

best_isozyme_for_smoment = argmax(
    smoment_isozyme_speed(gene_mass_function),
    my_isozyme_vector,
)
source
COBREXA.smoment_isozyme_speedMethod
smoment_isozyme_speed(isozyme::Isozyme, gene_product_molar_mass) -> Any

Compute a "score" for picking the most viable isozyme for make_smoment_model, based on maximum kcat divided by relative mass of the isozyme. This is used because sMOMENT algorithm can not handle multiple isozymes for one reaction.

source

Macro-generated functions and internal helpers

COBREXA.@_serialized_change_unwrapMacro
@_serialized_change_unwrap function

Creates a simple wrapper structure that calls the function transparently on the internal precached model. The internal type is returned (otherwise this would break the consistency of serialization).

source

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.

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

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

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

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

source
COBREXA.@_make_logging_tagMacro

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

source