index.html 63.5 KB
Newer Older
Documenter.jl's avatar
Documenter.jl committed
1
<!DOCTYPE html>
Documenter.jl's avatar
Documenter.jl committed
2
<html lang="en"><head><meta charset="UTF-8"/><meta name="viewport" content="width=device-width, initial-scale=1.0"/><title>Utilities · COBREXA.jl</title><link href="https://fonts.googleapis.com/css?family=Lato|Roboto+Mono" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/css/fontawesome.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/css/solid.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/css/brands.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/KaTeX/0.11.1/katex.min.css" rel="stylesheet" type="text/css"/><script>documenterBaseURL="../.."</script><script src="https://cdnjs.cloudflare.com/ajax/libs/require.js/2.3.6/require.min.js" data-main="../../assets/documenter.js"></script><script src="../../siteinfo.js"></script><script src="../../../versions.js"></script><link class="docs-theme-link" rel="stylesheet" type="text/css" href="../../assets/themes/documenter-dark.css" data-theme-name="documenter-dark" data-theme-primary-dark/><link class="docs-theme-link" rel="stylesheet" type="text/css" href="../../assets/themes/documenter-light.css" data-theme-name="documenter-light" data-theme-primary/><script src="../../assets/themeswap.js"></script><link href="../../assets/favicon.ico" rel="icon" type="image/x-icon"/></head><body><div id="documenter"><nav class="docs-sidebar"><a class="docs-logo" href="../../"><img class="docs-light-only" src="../../assets/logo.svg" alt="COBREXA.jl logo"/><img class="docs-dark-only" src="../../assets/logo-dark.svg" alt="COBREXA.jl logo"/></a><div class="docs-package-name"><span class="docs-autofit">COBREXA.jl</span></div><form class="docs-search" action="../../search/"><input class="docs-search-query" id="documenter-search-query" name="q" type="text" placeholder="Search docs"/></form><ul class="docs-menu"><li><a class="tocitem" href="../../">Home</a></li><li><span class="tocitem">User guide</span><ul><li><input class="collapse-toggle" id="menuitem-2-1" type="checkbox"/><label class="tocitem" for="menuitem-2-1"><span class="docs-label">Quickstart tutorials</span><i class="docs-chevron"></i></label><ul class="collapsed"><li><a class="tocitem" href="../../tutorials/">All tutorials</a></li><li><a class="tocitem" href="../../tutorials/1_loading/">Loading and converting model data</a></li><li><a class="tocitem" href="../../tutorials/2_analysis/">Basic analysis of constraint-based models</a></li><li><a class="tocitem" href="../../tutorials/3_hpc/">Distributed processing and HPC environments</a></li><li><a class="tocitem" href="../../tutorials/4_modifying/">Modifying and saving the models</a></li></ul></li><li><input class="collapse-toggle" id="menuitem-2-2" type="checkbox"/><label class="tocitem" for="menuitem-2-2"><span class="docs-label">Advanced tutorials</span><i class="docs-chevron"></i></label><ul class="collapsed"><li><a class="tocitem" href="../../advanced/">All advanced tutorials</a></li><li><a class="tocitem" href="../../advanced/1_variants/">Exploring many model variants</a></li><li><a class="tocitem" href="../../advanced/2_custom_model/">Working with custom models</a></li></ul></li><li><input class="collapse-toggle" id="menuitem-2-3" type="checkbox"/><label class="tocitem" for="menuitem-2-3"><span class="docs-label">Examples and notebooks</span><i class="docs-chevron"></i></label><ul class="collapsed"><li><a class="tocitem" href="../../notebooks/">All notebooks</a></li><li><a class="tocitem" href="../../notebooks/1_loading_converting_saving/">Loading, converting, and saving models</a></li><li><a class="tocitem" href="../../notebooks/2_finding_balance/">Finding balance and variability of constraint-based models</a></li><li><a class="tocitem" href="../../notebooks/3_basic_stdmodel_usage/">Basic usage of <code>StandardModel</code></a></li><li><a class="tocitem" href="../../notebooks/4_basic_core_coupled_usage/">Basic usage of <code>CoreModel</code> and <code>CoreModelCoupled</code></a></li><li><a class="tocitem" href="../../notebooks/5_basic_stdmodel_construction/">Model construction and modification</a></li><li><a class="tocitem" href="../../notebooks/6_screening/">Exploring model variants with <code>screen</code></a></li><li><a class="tocitem" href="../../notebooks/7_community_model/">Building and analysing a small community model</a></li><li><a class="tocitem" href="../../notebooks/8_custom_model/">Using a custom model data structure</a></li><li><a class="tocitem" href="../../notebooks/9_max_min_driving_force_analysis/">Maximum-minimum driving force analysis</a></li></ul></li></ul></li><li><span class="tocitem">Types and functions</span><ul><li><a class="tocitem" href="../">Contents</a></li><li><a class="tocitem" href="../analysis/">Analysis functions</a></li><li><a class="tocitem" href="../base/">Base functions</a></li><li><a class="tocitem" href="../io/">Input and output</a></li><li><a class="tocitem" href="../reconstruction/">Model construction functions</a></li><li><a class="tocitem" href="../types/">Types</a></li><li class="is-active"><a class="tocitem" href>Utilities</a><ul class="internal"><li><a class="tocitem" href="#Helper-functions"><span>Helper functions</span></a></li><li><a class="tocitem" href="#Macro-generated-functions-and-internal-helpers"><span>Macro-generated functions and internal helpers</span></a></li><li><a class="tocitem" href="#Logging-and-debugging-helpers"><span>Logging and debugging helpers</span></a></li></ul></li></ul></li><li><a class="tocitem" href="../../howToContribute/">How to contribute</a></li></ul><div class="docs-version-selector field has-addons"><div class="control"><span class="docs-label button is-static is-size-7">Version</span></div><div class="docs-selector control is-expanded"><div class="select is-fullwidth is-size-7"><select id="documenter-version-selector"></select></div></div></div></nav><div class="docs-main"><header class="docs-navbar"><nav class="breadcrumb"><ul class="is-hidden-mobile"><li><a class="is-disabled">Types and functions</a></li><li class="is-active"><a href>Utilities</a></li></ul><ul class="is-hidden-tablet"><li class="is-active"><a href>Utilities</a></li></ul></nav><div class="docs-right"><a class="docs-edit-link" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/master/docs/src/functions/utils.md" title="Edit on GitHub"><span class="docs-icon fab"></span><span class="docs-label is-hidden-touch">Edit on GitHub</span></a><a class="docs-settings-button fas fa-cog" id="documenter-settings-button" href="#" title="Settings"></a><a class="docs-sidebar-button fa fa-bars is-hidden-desktop" id="documenter-sidebar-button" href="#"></a></div></header><article class="content" id="documenter-page"><h1 id="Utilities"><a class="docs-heading-anchor" href="#Utilities">Utilities</a><a id="Utilities-1"></a><a class="docs-heading-anchor-permalink" href="#Utilities" title="Permalink"></a></h1><h2 id="Helper-functions"><a class="docs-heading-anchor" href="#Helper-functions">Helper functions</a><a id="Helper-functions-1"></a><a class="docs-heading-anchor-permalink" href="#Helper-functions" title="Permalink"></a></h2><article class="docstring"><header><a class="docstring-binding" id="COBREXA.ambiguously_identified_items-Tuple{Dict{String, Dict{String, Set{String}}}}" href="#COBREXA.ambiguously_identified_items-Tuple{Dict{String, Dict{String, Set{String}}}}"><code>COBREXA.ambiguously_identified_items</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">ambiguously_identified_items(
Documenter.jl's avatar
Documenter.jl committed
3
    index::Dict{String,Dict{String,[String]}},
Documenter.jl's avatar
Documenter.jl committed
4
)::Vector{String}</code></pre><p>Find items (genes, metabolites, ...) from the annotation index that are identified non-uniquely by at least one of their annotations.</p><p>This often indicates that the items are duplicate or miscategorized.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Annotation.jl#L44-L53">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.annotation_index-Tuple{AbstractDict{String}}" href="#COBREXA.annotation_index-Tuple{AbstractDict{String}}"><code>COBREXA.annotation_index</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">annotation_index(
Documenter.jl's avatar
Documenter.jl committed
5
6
    xs::AbstractDict{String};
    annotations = _annotations,
Documenter.jl's avatar
Documenter.jl committed
7
)::Dict{String,Dict{String,[String]}}</code></pre><p>Extract annotations from a dictionary of items <code>xs</code> and build an index that maps annotation &quot;kinds&quot; (e.g. <code>&quot;PubChem&quot;</code>) to the mapping from the annotations (e.g.  <code>&quot;COMPOUND_12345&quot;</code>) to item IDs that carry the annotations.</p><p>Function <code>annotations</code> is used to access the <code>Annotations</code> object in the dictionary values.</p><p>This is extremely useful for finding items by annotation data.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Annotation.jl#L6-L20">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.check_duplicate_reaction-Tuple{Reaction, OrderedCollections.OrderedDict{String, Reaction}}" href="#COBREXA.check_duplicate_reaction-Tuple{Reaction, OrderedCollections.OrderedDict{String, Reaction}}"><code>COBREXA.check_duplicate_reaction</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">check_duplicate_reaction(rxn::Reaction, rxns::Dict{String, Reaction}; only_metabolites=true)</code></pre><p>Check if <code>rxn</code> already exists in <code>rxns</code> but has another <code>id</code>. If <code>only_metabolites</code> is <code>true</code> then only the metabolite <code>id</code>s are checked. Otherwise, compares metabolite <code>id</code>s and the absolute value of their stoichiometric coefficients to those of <code>rxn</code>. If <code>rxn</code> has the same reaction equation as another reaction in <code>rxns</code>, the return the <code>id</code>. Otherwise return <code>nothing</code>.</p><p>See also: <a href="#COBREXA.reaction_mass_balanced-Tuple{StandardModel, String}"><code>reaction_mass_balanced</code></a></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Reaction.jl#L1-L11">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.is_boundary-Tuple{Dict{String, Float64}}" href="#COBREXA.is_boundary-Tuple{Dict{String, Float64}}"><code>COBREXA.is_boundary</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">is_boundary(rxn_dict::Dict{String, Float64})</code></pre><p>Return true if the reaction denoted by <code>rxn_dict</code> is a boundary reaction, otherwise return false. Checks if on boundary by inspecting the number of metabolites in <code>rxn_dict</code>. Boundary reactions have only one metabolite, e.g. an exchange reaction, or a sink/demand reaction.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Reaction.jl#L41-L47">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.reaction_atom_balance-Tuple{StandardModel, Dict{String, Float64}}" href="#COBREXA.reaction_atom_balance-Tuple{StandardModel, Dict{String, Float64}}"><code>COBREXA.reaction_atom_balance</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">reaction_atom_balance(model::StandardModel, rxn)</code></pre><p>Returns a dictionary mapping the stoichiometry of atoms through a single reaction. Uses the metabolite information in <code>model</code> to determine the mass balance. Accepts a reaction dictionary, a reaction string id or a <code>Reaction</code> as an argument for <code>rxn</code>.</p><p>See also: <a href="#COBREXA.reaction_mass_balanced-Tuple{StandardModel, String}"><code>reaction_mass_balanced</code></a></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Reaction.jl#L59-L67">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.reaction_mass_balanced-Tuple{StandardModel, String}" href="#COBREXA.reaction_mass_balanced-Tuple{StandardModel, String}"><code>COBREXA.reaction_mass_balanced</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">reaction_mass_balanced(model::StandardModel, rxn)</code></pre><p>Checks if <code>rxn</code> is atom balanced. Returns a boolean for whether the reaction is balanced, and the associated balance of atoms for convenience (useful if not balanced). Calls <code>reaction_atom_balance</code> internally.</p><p>See also: <a href="#COBREXA.check_duplicate_reaction-Tuple{Reaction, OrderedCollections.OrderedDict{String, Reaction}}"><code>check_duplicate_reaction</code></a>, <a href="#COBREXA.reaction_atom_balance-Tuple{StandardModel, Dict{String, Float64}}"><code>reaction_atom_balance</code></a></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Reaction.jl#L89-L97">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.stoichiometry_string-Tuple{Any}" href="#COBREXA.stoichiometry_string-Tuple{Any}"><code>COBREXA.stoichiometry_string</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">stoichiometry_string(rxn_dict::Dict{String, Float64}; format_id = x -&gt; x)</code></pre><p>Return the reaction equation as a string. The metabolite strings can be manipulated by setting <code>format_id</code>.</p><p><strong>Example</strong></p><pre><code class="language-none">julia&gt; req = Dict(&quot;coa_c&quot; =&gt; -1, &quot;for_c&quot; =&gt; 1, &quot;accoa_c&quot; =&gt; 1, &quot;pyr_c&quot; =&gt; -1)
Documenter.jl's avatar
Documenter.jl committed
8
9
10
11
julia&gt; stoichiometry_string(req)
&quot;coa_c + pyr_c = for_c + accoa_c&quot;

julia&gt; stoichiometry_string(req; format_id = x -&gt; x[1:end-2])
Documenter.jl's avatar
Documenter.jl committed
12
&quot;coa + pyr = for + accoa&quot;</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Reaction.jl#L107-L122">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.stoichiometry_string-Tuple{Reaction}" href="#COBREXA.stoichiometry_string-Tuple{Reaction}"><code>COBREXA.stoichiometry_string</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">stoichiometry_string(rxn::Reaction; kwargs)</code></pre><p>Alternative of <a href="#COBREXA.stoichiometry_string-Tuple{Any}"><code>stoichiometry_string</code></a> take takes a <code>Reaction</code> as an argument.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Reaction.jl#L132-L136">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.serialize_model-Tuple{Serialized, String}" href="#COBREXA.serialize_model-Tuple{Serialized, String}"><code>COBREXA.serialize_model</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">serialize_model(model::Serialized, filename::String)::Serialized</code></pre><p>Specialization of <a href="#COBREXA.serialize_model-Tuple{Serialized, String}"><code>serialize_model</code></a> that prevents nested serialization of already-serialized models.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Serialized.jl#L24-L29">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.serialize_model-Union{Tuple{MM}, Tuple{MM, String}} where MM&lt;:MetabolicModel" href="#COBREXA.serialize_model-Union{Tuple{MM}, Tuple{MM, String}} where MM&lt;:MetabolicModel"><code>COBREXA.serialize_model</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">serialize_model(model::MM, filename::String)::Serialized{MM} where {MM&lt;:MetabolicModel}</code></pre><p>Serialize the <code>model</code> to file <code>filename</code>, returning a <a href="../types/#COBREXA.Serialized"><code>Serialized</code></a> model that can be loaded back transparently by <a href="../types/#COBREXA.precache!-Tuple{MetabolicModel}"><code>precache!</code></a>. The result does <em>not</em> 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.</p><p>The benefit of using this over &quot;raw&quot; <code>Serialization.serialize</code> is that the resulting <code>Serialized</code> model will reload itself automatically with <a href="../types/#COBREXA.precache!-Tuple{MetabolicModel}"><code>precache!</code></a> at first usage, which needs to be done manually when using the <code>Serialization</code> package directly.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/Serialized.jl#L2-L15">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.copy-Tuple{Gene}" href="#Base.copy-Tuple{Gene}"><code>Base.copy</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">Base.copy(g::Gene)</code></pre><p>Shallow copy of a <a href="../types/#COBREXA.Gene"><code>Gene</code></a></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/StandardModel.jl#L45-L49">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.copy-Tuple{Metabolite}" href="#Base.copy-Tuple{Metabolite}"><code>Base.copy</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">Base.copy(m::Metabolite)</code></pre><p>Shallow copy of a <a href="../types/#COBREXA.Metabolite"><code>Metabolite</code></a></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/StandardModel.jl#L31-L35">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.copy-Tuple{Reaction}" href="#Base.copy-Tuple{Reaction}"><code>Base.copy</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">Base.copy(r::Reaction)</code></pre><p>Shallow copy of a <a href="../types/#COBREXA.Reaction"><code>Reaction</code></a></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/StandardModel.jl#L14-L18">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.copy-Tuple{StandardModel}" href="#Base.copy-Tuple{StandardModel}"><code>Base.copy</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">Base.copy(m::StandardModel)</code></pre><p>Shallow copy of a <a href="../types/#COBREXA.StandardModel"><code>StandardModel</code></a></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/StandardModel.jl#L2-L6">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gamma_bounds-Tuple{Any}" href="#COBREXA.gamma_bounds-Tuple{Any}"><code>COBREXA.gamma_bounds</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gamma_bounds(gamma)</code></pre><p>A bounds-generating function for <a href="../analysis/#COBREXA.flux_variability_analysis-Tuple{MetabolicModel, Any}"><code>flux_variability_analysis</code></a> that limits the objective value to be at least <code>gamma*Z₀</code>, as usual in COBRA packages. Use as the <code>bounds</code> argument:</p><pre><code class="language-none">flux_variability_analysis(model, some_optimizer; bounds = gamma_bounds(0.9))</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/bounds.jl#L1-L10">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.objective_bounds-Tuple{Any}" href="#COBREXA.objective_bounds-Tuple{Any}"><code>COBREXA.objective_bounds</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">objective_bounds(tolerance)</code></pre><p>A bounds-generating function for <a href="../analysis/#COBREXA.flux_variability_analysis-Tuple{MetabolicModel, Any}"><code>flux_variability_analysis</code></a> that limits the objective value to a small multiple of Z₀. Use as <code>bounds</code> argument, similarly to <a href="#COBREXA.gamma_bounds-Tuple{Any}"><code>gamma_bounds</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/bounds.jl#L13-L19">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._parse_formula-Tuple{String}" href="#COBREXA._parse_formula-Tuple{String}"><code>COBREXA._parse_formula</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_parse_formula(f::String)::MetaboliteFormula</code></pre><p>Parse a formula in format <code>C2H6O</code> into a <a href="../types/#COBREXA.MetaboliteFormula"><code>MetaboliteFormula</code></a>, which is basically a dictionary of atom counts in the molecule.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/chemical_formulas.jl#L2-L7">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._unparse_formula-Tuple{Dict{String, Int64}}" href="#COBREXA._unparse_formula-Tuple{Dict{String, Int64}}"><code>COBREXA._unparse_formula</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_unparse_formula(f::MetaboliteFormula)::String</code></pre><p>Format <a href="../types/#COBREXA.MetaboliteFormula"><code>MetaboliteFormula</code></a> to <code>String</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/chemical_formulas.jl#L19-L23">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gene_product_dict-Tuple{GeckoModel, Any}" href="#COBREXA.gene_product_dict-Tuple{GeckoModel, Any}"><code>COBREXA.gene_product_dict</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gene_product_dict(model::GeckoModel, opt_model)</code></pre><p>Return a dictionary mapping protein molar concentrations to their ids. The argument <code>opt_model</code> is a solved optimization problem, typically returned by <a href="../analysis/#COBREXA.flux_balance_analysis-Union{Tuple{M}, Tuple{M, Any}} where M&lt;:MetabolicModel"><code>flux_balance_analysis</code></a>. See <a href="../base/#COBREXA.flux_dict-Tuple{MetabolicModel, Any}"><code>flux_dict</code></a> for the corresponding function that returns a dictionary of solved fluxes.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/enzymes.jl#L1-L8">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gene_product_dict-Tuple{GeckoModel}" href="#COBREXA.gene_product_dict-Tuple{GeckoModel}"><code>COBREXA.gene_product_dict</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gene_product_dict(model::GeckoModel)</code></pre><p>A pipe-able variant of <a href="#COBREXA.gene_product_dict-Tuple{GeckoModel, Any}"><code>gene_product_dict</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/enzymes.jl#L13-L17">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gene_product_mass-Tuple{SMomentModel, Any}" href="#COBREXA.gene_product_mass-Tuple{SMomentModel, Any}"><code>COBREXA.gene_product_mass</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gene_product_mass(model::SMomentModel)</code></pre><p>Extract the total mass utilization in a solved <a href="../types/#COBREXA.SMomentModel"><code>SMomentModel</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/enzymes.jl#L42-L46">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gene_product_mass-Tuple{SMomentModel}" href="#COBREXA.gene_product_mass-Tuple{SMomentModel}"><code>COBREXA.gene_product_mass</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gene_product_mass(model::SMomentModel)</code></pre><p>A pipe-able variant of <a href="#COBREXA.gene_product_mass-Tuple{SMomentModel, Any}"><code>gene_product_mass</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/enzymes.jl#L51-L56">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gene_product_mass_group_dict-Tuple{GeckoModel, Any}" href="#COBREXA.gene_product_mass_group_dict-Tuple{GeckoModel, Any}"><code>COBREXA.gene_product_mass_group_dict</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gene_product_mass_group_dict(model::GeckoModel, opt_model)</code></pre><p>Extract the mass utilization in mass groups from a solved <a href="../types/#COBREXA.GeckoModel"><code>GeckoModel</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/enzymes.jl#L20-L24">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gene_product_mass_group_dict-Tuple{GeckoModel}" href="#COBREXA.gene_product_mass_group_dict-Tuple{GeckoModel}"><code>COBREXA.gene_product_mass_group_dict</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gene_product_mass_group_dict(model::GeckoModel)</code></pre><p>A pipe-able variant of <a href="#COBREXA.gene_product_mass_group_dict-Tuple{GeckoModel, Any}"><code>gene_product_mass_group_dict</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/enzymes.jl#L34-L38">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.atom_fluxes-Tuple{MetabolicModel, Dict{String, Float64}}" href="#COBREXA.atom_fluxes-Tuple{MetabolicModel, Dict{String, Float64}}"><code>COBREXA.atom_fluxes</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">atom_fluxes(model::MetabolicModel, reaction_fluxes::Dict{String, Float64})</code></pre><p>Return a dictionary mapping the flux of atoms across a flux solution given by <code>reaction_fluxes</code> using the reactions in <code>model</code> to determine the appropriate stoichiometry.</p><p>Note, this function ignores metabolites with no formula assigned to them, no error message will be displayed.</p><p>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.</p><p><strong>Example</strong></p><pre><code class="language-none"># Find flux of Carbon through all metabolic reactions except the biomass reaction
Documenter.jl's avatar
Documenter.jl committed
13
delete!(fluxes, &quot;BIOMASS_Ecoli_core_w_GAM&quot;)
Documenter.jl's avatar
Documenter.jl committed
14
atom_fluxes(model, fluxes)[&quot;C&quot;]</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/fluxes.jl#L35-L53">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.metabolite_fluxes-Tuple{MetabolicModel, Dict{String, Float64}}" href="#COBREXA.metabolite_fluxes-Tuple{MetabolicModel, Dict{String, Float64}}"><code>COBREXA.metabolite_fluxes</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">metabolite_fluxes(model::MetabolicModel, flux_dict::Dict{String, Float64})</code></pre><p>Return two dictionaries of metabolite <code>id</code>s mapped to reactions that consume or produce them, given the flux distribution supplied in <code>flux_dict</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/fluxes.jl#L1-L6">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._gecko_gene_product_coupling-Tuple{GeckoModel}" href="#COBREXA._gecko_gene_product_coupling-Tuple{GeckoModel}"><code>COBREXA._gecko_gene_product_coupling</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_gecko_gene_product_coupling(model::GeckoModel)</code></pre><p>Compute the part of the coupling for GeckoModel that limits the amount of each kind of protein available.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gecko.jl#L54-L59">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._gecko_mass_group_coupling-Tuple{GeckoModel}" href="#COBREXA._gecko_mass_group_coupling-Tuple{GeckoModel}"><code>COBREXA._gecko_mass_group_coupling</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_gecko_mass_group_coupling(model::GeckoModel)</code></pre><p>Compute the part of the coupling for <a href="../types/#COBREXA.GeckoModel"><code>GeckoModel</code></a> that limits the total mass of each group of gene products.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gecko.jl#L75-L80">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._gecko_reaction_column_reactions-Tuple{Any, Any}" href="#COBREXA._gecko_reaction_column_reactions-Tuple{Any, Any}"><code>COBREXA._gecko_reaction_column_reactions</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_gecko_reaction_column_reactions(columns, inner)</code></pre><p>Helper method that doesn&#39;t require the whole <a href="../types/#COBREXA.GeckoModel"><code>GeckoModel</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gecko.jl#L21-L25">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._gecko_reaction_column_reactions-Tuple{GeckoModel}" href="#COBREXA._gecko_reaction_column_reactions-Tuple{GeckoModel}"><code>COBREXA._gecko_reaction_column_reactions</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_gecko_reaction_column_reactions(model::GeckoModel)</code></pre><p>Retrieve a utility mapping between reactions and split reactions; rows correspond to &quot;original&quot; reactions, columns correspond to &quot;split&quot; reactions.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gecko.jl#L12-L17">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._gecko_reaction_coupling-Tuple{GeckoModel}" href="#COBREXA._gecko_reaction_coupling-Tuple{GeckoModel}"><code>COBREXA._gecko_reaction_coupling</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_gecko_reaction_coupling(model::GeckoModel)</code></pre><p>Compute the part of the coupling for <a href="../types/#COBREXA.GeckoModel"><code>GeckoModel</code></a> that limits the &quot;arm&quot; reactions (which group the individual split unidirectional reactions).</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gecko.jl#L34-L39">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._gecko_reaction_name-Tuple{String, Int64, Int64}" href="#COBREXA._gecko_reaction_name-Tuple{String, Int64, Int64}"><code>COBREXA._gecko_reaction_name</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_gecko_reaction_name(original_name::String, direction::Int)</code></pre><p>Internal helper for systematically naming reactions in <a href="../types/#COBREXA.GeckoModel"><code>GeckoModel</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gecko.jl#L2-L6">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._parse_grr-Tuple{SBML.GeneProductAssociation}" href="#COBREXA._parse_grr-Tuple{SBML.GeneProductAssociation}"><code>COBREXA._parse_grr</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_parse_grr(gpa::SBML.GeneProductAssociation)::GeneAssociation</code></pre><p>Parse <code>SBML.GeneProductAssociation</code> structure to the simpler GeneAssociation. The input must be (implicitly) in a positive DNF.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gene_associations.jl#L2-L7">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._parse_grr-Tuple{String}" href="#COBREXA._parse_grr-Tuple{String}"><code>COBREXA._parse_grr</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_parse_grr(s::String)::GeneAssociation</code></pre><p>Parse a DNF gene association rule in format <code>(YIL010W and YLR043C) or (YIL010W and YGR209C)</code> to <code>GeneAssociation. Also accepts</code>OR<code>,</code>|<code>,</code>||<code>,</code>AND<code>,</code>&amp;<code>, and</code>&amp;&amp;`.</p><p><strong>Example</strong></p><pre><code class="language-none">julia&gt; _parse_grr(&quot;(YIL010W and YLR043C) or (YIL010W and YGR209C)&quot;)
Documenter.jl's avatar
Documenter.jl committed
15
16
2-element Array{Array{String,1},1}:
 [&quot;YIL010W&quot;, &quot;YLR043C&quot;]
Documenter.jl's avatar
Documenter.jl committed
17
 [&quot;YIL010W&quot;, &quot;YGR209C&quot;]</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gene_associations.jl#L37-L51">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._parse_grr_to_sbml-Tuple{String}" href="#COBREXA._parse_grr_to_sbml-Tuple{String}"><code>COBREXA._parse_grr_to_sbml</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_parse_grr_to_sbml(str::String)::Maybe{SBML.GeneProductAssociation}</code></pre><p>Internal helper for parsing the string GRRs into SBML data structures. More general than <a href="#COBREXA._parse_grr-Tuple{SBML.GeneProductAssociation}"><code>_parse_grr</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gene_associations.jl#L54-L59">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._unparse_grr-Tuple{Type{SBML.GeneProductAssociation}, Vector{Vector{String}}}" href="#COBREXA._unparse_grr-Tuple{Type{SBML.GeneProductAssociation}, Vector{Vector{String}}}"><code>COBREXA._unparse_grr</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_unparse_grr(
Documenter.jl's avatar
Documenter.jl committed
18
19
    ::Type{SBML.GeneProductAssociation},
    x::GeneAssociation,
Documenter.jl's avatar
Documenter.jl committed
20
21
)::SBML.GeneAssociation</code></pre><p>Convert a GeneAssociation to the corresponding <code>SBML.jl</code> structure.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gene_associations.jl#L22-L29">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._unparse_grr-Tuple{Type{String}, Vector{Vector{String}}}" href="#COBREXA._unparse_grr-Tuple{Type{String}, Vector{Vector{String}}}"><code>COBREXA._unparse_grr</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">unparse_grr(grr::Vector{Vector{Gene}}</code></pre><p>Converts a nested string gene reaction array  back into a gene reaction rule string.</p><p><strong>Example</strong></p><pre><code class="language-none">julia&gt; _unparse_grr(String, [[&quot;YIL010W&quot;, &quot;YLR043C&quot;], [&quot;YIL010W&quot;, &quot;YGR209C&quot;]])
&quot;(YIL010W and YLR043C) or (YIL010W and YGR209C)&quot;</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/gene_associations.jl#L119-L130">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._guesskey-Tuple{Any, Any}" href="#COBREXA._guesskey-Tuple{Any, Any}"><code>COBREXA._guesskey</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_guesskey(ks, possibilities)</code></pre><p>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 <code>src/base/constants.jl</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/guesskey.jl#L2-L9">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.gets-Tuple{Any, Any, Any}" href="#COBREXA.gets-Tuple{Any, Any, Any}"><code>COBREXA.gets</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">gets(collection, fail, keys)</code></pre><p>Return <code>fail</code> if key in <code>keys</code> is not in <code>collection</code>, otherwise return <code>collection[key]</code>. Useful if may different keys need to be tried due to non-standardized model formats.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/guesskey.jl#L24-L30">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.find_biomass_reaction_ids-Tuple{MetabolicModel}" href="#COBREXA.find_biomass_reaction_ids-Tuple{MetabolicModel}"><code>COBREXA.find_biomass_reaction_ids</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">find_biomass_reaction_ids(m::MetabolicModel; kwargs...)</code></pre><p>Shortcut for finding biomass reaction identifiers in a model; arguments are forwarded to <a href="#COBREXA.looks_like_biomass_reaction-Tuple{String}"><code>looks_like_biomass_reaction</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L92-L97">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.find_biomass_reactions-Tuple{MetabolicModel}" href="#COBREXA.find_biomass_reactions-Tuple{MetabolicModel}"><code>COBREXA.find_biomass_reactions</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">find_biomass_reactions(m::MetabolicModel; kwargs...)</code></pre><p>Shortcut for finding biomass reaction indexes in a model; arguments are forwarded to <a href="#COBREXA.looks_like_biomass_reaction-Tuple{String}"><code>looks_like_biomass_reaction</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L83-L88">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.find_exchange_reaction_ids-Tuple{MetabolicModel}" href="#COBREXA.find_exchange_reaction_ids-Tuple{MetabolicModel}"><code>COBREXA.find_exchange_reaction_ids</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">find_exchange_reaction_ids(m::MetabolicModel; kwargs...)</code></pre><p>Shortcut for finding exchange reaction identifiers in a model; arguments are forwarded to <a href="#COBREXA.looks_like_exchange_reaction-Tuple{String}"><code>looks_like_exchange_reaction</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L46-L51">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.find_exchange_reactions-Tuple{MetabolicModel}" href="#COBREXA.find_exchange_reactions-Tuple{MetabolicModel}"><code>COBREXA.find_exchange_reactions</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">find_exchange_reactions(m::MetabolicModel; kwargs...)</code></pre><p>Shortcut for finding exchange reaction indexes in a model; arguments are forwarded to <a href="#COBREXA.looks_like_exchange_reaction-Tuple{String}"><code>looks_like_exchange_reaction</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L37-L42">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.find_extracellular_metabolite_ids-Tuple{MetabolicModel}" href="#COBREXA.find_extracellular_metabolite_ids-Tuple{MetabolicModel}"><code>COBREXA.find_extracellular_metabolite_ids</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">find_extracellular_metabolite_ids(m::MetabolicModel; kwargs...)</code></pre><p>Shortcut for finding extracellular metabolite identifiers in a model; arguments are forwarded to <a href="#COBREXA.looks_like_extracellular_metabolite-Tuple{String}"><code>looks_like_extracellular_metabolite</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L131-L136">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.find_extracellular_metabolites-Tuple{MetabolicModel}" href="#COBREXA.find_extracellular_metabolites-Tuple{MetabolicModel}"><code>COBREXA.find_extracellular_metabolites</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">find_extracellular_metabolites(m::MetabolicModel; kwargs...)</code></pre><p>Shortcut for finding extracellular metabolite indexes in a model; arguments are forwarded to <a href="#COBREXA.looks_like_extracellular_metabolite-Tuple{String}"><code>looks_like_extracellular_metabolite</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L123-L128">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.looks_like_biomass_reaction-Tuple{String}" href="#COBREXA.looks_like_biomass_reaction-Tuple{String}"><code>COBREXA.looks_like_biomass_reaction</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">looks_like_biomass_reaction(rxn_id::String;
Documenter.jl's avatar
Documenter.jl committed
22
23
24
25
    exclude_exchanges = false,
    exchange_prefixes = _constants.exchange_prefixes,
    biomass_strings = _constants.biomass_strings,
)::Bool</code></pre><p>A predicate that matches reaction identifiers that look like biomass reactions. Biomass reactions are identified by looking for occurences of <code>biomass_strings</code> in the reaction id. If <code>exclude_exchanges</code> is set, the strings that look like exchanges (from <a href="#COBREXA.looks_like_exchange_reaction-Tuple{String}"><code>looks_like_exchange_reaction</code></a>) will not match.</p><p><strong>Example</strong></p><pre><code class="language-none">filter(looks_like_biomass_reaction, reactions(model)) # returns strings
Documenter.jl's avatar
Documenter.jl committed
26
findall(looks_like_biomass_reaction, reactions(model)) # returns indices</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L55-L72">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.looks_like_exchange_reaction-Tuple{String}" href="#COBREXA.looks_like_exchange_reaction-Tuple{String}"><code>COBREXA.looks_like_exchange_reaction</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">looks_like_exchange_reaction(rxn_id::String;
Documenter.jl's avatar
Documenter.jl committed
27
28
29
30
31
32
33
34
35
    exclude_biomass = false,
    biomass_strings = _constants.biomass_strings,
    exchange_prefixes = _constants.exchange_prefixes,
)</code></pre><p>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 <code>exchange_prefixes</code> and biomass reactions are identified by looking for occurences of <code>biomass_strings</code> in the reaction id.</p><p>Also see <a href="#COBREXA.find_exchange_reactions-Tuple{MetabolicModel}"><code>find_exchange_reactions</code></a>.</p><p><strong>Example</strong></p><pre><code class="language-none">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&#39;s arguments
# using an anonymous function
findall(x -&gt; looks_like_exchange_reaction(x; exclude_biomass=true), reactions(model)) # returns indices
Documenter.jl's avatar
Documenter.jl committed
36
filter(x -&gt; looks_like_exchange_reaction(x; exclude_biomass=true), reactions(model)) # returns Strings</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L1-L26">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.looks_like_extracellular_metabolite-Tuple{String}" href="#COBREXA.looks_like_extracellular_metabolite-Tuple{String}"><code>COBREXA.looks_like_extracellular_metabolite</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">looks_like_extracellular_metabolite(rxn_id::String;
Documenter.jl's avatar
Documenter.jl committed
37
    extracellular_suffixes = _constants.extracellular_suffixes,
Documenter.jl's avatar
Documenter.jl committed
38
)::Bool</code></pre><p>A predicate that matches metabolite identifiers that look like they are extracellular metabolites. Extracellular metabolites are identified by <code>extracellular_suffixes</code> at the end of the metabolite id.</p><p><strong>Example</strong></p><pre><code class="language-none">filter(looks_like_extracellular_metabolite, metabolites(model)) # returns strings
Documenter.jl's avatar
Documenter.jl committed
39
findall(looks_like_extracellular_metabolite, metabolites(model)) # returns indices</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/looks_like.jl#L101-L115">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._smoment_column_reactions-Tuple{SMomentModel}" href="#COBREXA._smoment_column_reactions-Tuple{SMomentModel}"><code>COBREXA._smoment_column_reactions</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_smoment_column_reactions(model::SMomentModel)</code></pre><p>Retrieve a utility mapping between reactions and split reactions; rows correspond to &quot;original&quot; reactions, columns correspond to &quot;split&quot; reactions.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/smoment.jl#L11-L16">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._smoment_reaction_name-Tuple{String, Int64}" href="#COBREXA._smoment_reaction_name-Tuple{String, Int64}"><code>COBREXA._smoment_reaction_name</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_smoment_reaction_name(original_name::String, direction::Int)</code></pre><p>Internal helper for systematically naming reactions in <a href="../types/#COBREXA.SMomentModel"><code>SMomentModel</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/smoment.jl#L2-L6">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.smoment_isozyme_speed-Tuple{Function}" href="#COBREXA.smoment_isozyme_speed-Tuple{Function}"><code>COBREXA.smoment_isozyme_speed</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">smoment_isozyme_speed(gene_product_molar_mass::Function)</code></pre><p>A piping- and argmax-friendly overload of <a href="#COBREXA.smoment_isozyme_speed-Tuple{Function}"><code>smoment_isozyme_speed</code></a>.</p><p><strong>Example</strong></p><pre><code class="language-none">gene_mass_function = gid -&gt; 1.234
Documenter.jl's avatar
Documenter.jl committed
40
41
42
43

best_isozyme_for_smoment = argmax(
    smoment_isozyme_speed(gene_mass_function),
    my_isozyme_vector,
Documenter.jl's avatar
Documenter.jl committed
44
45
)</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/smoment.jl#L39-L53">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.smoment_isozyme_speed-Tuple{Isozyme, Any}" href="#COBREXA.smoment_isozyme_speed-Tuple{Isozyme, Any}"><code>COBREXA.smoment_isozyme_speed</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">smoment_isozyme_speed(isozyme::Isozyme, gene_product_molar_mass)</code></pre><p>Compute a &quot;score&quot; for picking the most viable isozyme for <a href="../analysis/#COBREXA.make_smoment_model-Tuple{MetabolicModel}"><code>make_smoment_model</code></a>, 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.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/utils/smoment.jl#L25-L32">source</a></section></article><h2 id="Macro-generated-functions-and-internal-helpers"><a class="docs-heading-anchor" href="#Macro-generated-functions-and-internal-helpers">Macro-generated functions and internal helpers</a><a id="Macro-generated-functions-and-internal-helpers-1"></a><a class="docs-heading-anchor-permalink" href="#Macro-generated-functions-and-internal-helpers" title="Permalink"></a></h2><article class="docstring"><header><a class="docstring-binding" id="COBREXA.@_change_bounds_fn-Tuple{Any, Any, Vararg{Any}}" href="#COBREXA.@_change_bounds_fn-Tuple{Any, Any, Vararg{Any}}"><code>COBREXA.@_change_bounds_fn</code></a><span class="docstring-category">Macro</span></header><section><div><pre><code class="language-julia">@_change_bounds_fn ModelType IdxType [plural] [inplace] begin ... end</code></pre><p>A helper for creating simple bounds-changing function similar to <a href="../reconstruction/#COBREXA.change_bounds-Tuple{CoreModel, AbstractVector{Int64}}"><code>change_bounds</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/macros/change_bounds.jl#L1-L6">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA._inherit_model_methods_impl-Tuple{Any, Symbol, Any, Any, Any, Vararg{Any}}" href="#COBREXA._inherit_model_methods_impl-Tuple{Any, Symbol, Any, Any, Any, Vararg{Any}}"><code>COBREXA._inherit_model_methods_impl</code></a><span class="docstring-category">Method</span></header><section><div><pre><code class="language-julia">_inherit_model_methods_impl(mtype::Symbol, arglist, access, fwdlist, fns...)</code></pre><p>A helper backend for <a href="#COBREXA.@_inherit_model_methods-Tuple{Symbol, Any, Symbol, Any, Vararg{Any}}"><code>@_inherit_model_methods</code></a> and <a href="#COBREXA.@_inherit_model_methods_fn-Tuple{Symbol, Any, Any, Any, Vararg{Any}}"><code>@_inherit_model_methods_fn</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/macros/model_wrapper.jl#L2-L7">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.@_inherit_model_methods-Tuple{Symbol, Any, Symbol, Any, Vararg{Any}}" href="#COBREXA.@_inherit_model_methods-Tuple{Symbol, Any, Symbol, Any, Vararg{Any}}"><code>COBREXA.@_inherit_model_methods</code></a><span class="docstring-category">Macro</span></header><section><div><pre><code class="language-julia">@_inherit_model_methods</code></pre><p>Generates trivial accessor functions listed in <code>fns</code> for a model that is wrapped in type <code>mtype</code> as field <code>member</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/macros/model_wrapper.jl#L40-L45">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.@_inherit_model_methods_fn-Tuple{Symbol, Any, Any, Any, Vararg{Any}}" href="#COBREXA.@_inherit_model_methods_fn-Tuple{Symbol, Any, Any, Any, Vararg{Any}}"><code>COBREXA.@_inherit_model_methods_fn</code></a><span class="docstring-category">Macro</span></header><section><div><pre><code class="language-julia">@_inherit_model_methods_fn</code></pre><p>A more generic version of <a href="#COBREXA.@_inherit_model_methods-Tuple{Symbol, Any, Symbol, Any, Vararg{Any}}"><code>@_inherit_model_methods</code></a> that accesses the &quot;inner&quot; model using an accessor function name.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/macros/model_wrapper.jl#L57-L62">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.@_remove_fn-Tuple{Any, Any, Any, Vararg{Any}}" href="#COBREXA.@_remove_fn-Tuple{Any, Any, Any, Vararg{Any}}"><code>COBREXA.@_remove_fn</code></a><span class="docstring-category">Macro</span></header><section><div><pre><code class="language-julia">@ _remove_fn objname ModelType IndexType [plural] [inplace] begin ... end</code></pre><p>A helper for creating functions that follow the <code>remove_objname</code> template, such as <a href="../reconstruction/#COBREXA.remove_metabolites-Tuple{CoreModel, AbstractVector{Int64}}"><code>remove_metabolites</code></a> and <a href="../reconstruction/#COBREXA.remove_reaction-Tuple{CoreModel, Int64}"><code>remove_reaction</code></a>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/macros/remove_item.jl#L2-L7">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.@_serialized_change_unwrap-Tuple{Symbol}" href="#COBREXA.@_serialized_change_unwrap-Tuple{Symbol}"><code>COBREXA.@_serialized_change_unwrap</code></a><span class="docstring-category">Macro</span></header><section><div><pre><code class="language-julia">@_serialized_change_unwrap function</code></pre><p>Creates a simple wrapper structure that calls the <code>function</code> transparently on the internal precached model. The internal type is returned (otherwise this would break the consistency of serialization).</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/macros/serialized.jl#L1-L7">source</a></section></article><h2 id="Logging-and-debugging-helpers"><a class="docs-heading-anchor" href="#Logging-and-debugging-helpers">Logging and debugging helpers</a><a id="Logging-and-debugging-helpers-1"></a><a class="docs-heading-anchor-permalink" href="#Logging-and-debugging-helpers" title="Permalink"></a></h2><article class="docstring"><header><a class="docstring-binding" id="COBREXA.log_io" href="#COBREXA.log_io"><code>COBREXA.log_io</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">log_io(enable::Bool=true)</code></pre><p>Enable (default) or disable (by passing <code>false</code>) output of messages and warnings from model input/output.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/logging/log.jl#L40-L44">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.log_models" href="#COBREXA.log_models"><code>COBREXA.log_models</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">log_models(enable::Bool=true)</code></pre><p>Enable (default) or disable (by passing <code>false</code>) output of model-related messages.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/logging/log.jl#L40-L44">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.log_perf" href="#COBREXA.log_perf"><code>COBREXA.log_perf</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">log_perf(enable::Bool=true)</code></pre><p>Enable (default) or disable (by passing <code>false</code>) output of performance-related tracing information.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/logging/log.jl#L40-L44">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="COBREXA.@_make_logging_tag-Tuple{Symbol, String}" href="#COBREXA.@_make_logging_tag-Tuple{Symbol, String}"><code>COBREXA.@_make_logging_tag</code></a><span class="docstring-category">Macro</span></header><section><div><pre><code class="language-julia">macro _make_logging_group(sym::Symbol, doc::String)</code></pre><p>This creates a group of functions that allow masking out topic-related logging actions. A call that goes as follows:</p><pre><code class="language-none">@_make_logging_tag XYZ</code></pre><p>creates the following tools:</p><ul><li>global variable <code>_XYZ_log_enabled</code> defaulted to false</li><li>function <code>log_XYZ</code> that can be called to turn the logging on/off</li><li>a masking macro <code>@_XYZ_log</code> that can be prepended to commands that should only happen if the logging of tag XYZ is enabled.</li></ul><p>The masking macro is then used as follows:</p><pre><code class="language-none">@_XYZ_log @info &quot;This is the extra verbose information you wanted!&quot; a b c</code></pre><p>The user can direct logging with these:</p><pre><code class="language-none">log_XYZ()
log_XYZ(false)</code></pre><p><code>doc</code> 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 <code>&quot;X, Y and Z-related messages&quot;</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/lcsb-biocore/COBREXA.jl/blob/77603c3b626a5edd84a1fbd03d11477955788b00/src/base/logging/log.jl#L2-L29">source</a></section></article></article><nav class="docs-footer"><a class="docs-footer-prevpage" href="../types/">« Types</a><a class="docs-footer-nextpage" href="../../howToContribute/">How to contribute »</a><div class="flexbox-break"></div><p class="footer-message">Powered by <a href="https://github.com/JuliaDocs/Documenter.jl">Documenter.jl</a> and the <a href="https://julialang.org/">Julia Programming Language</a>.</p></nav></div><div class="modal" id="documenter-settings"><div class="modal-background"></div><div class="modal-card"><header class="modal-card-head"><p class="modal-card-title">Settings</p><button class="delete"></button></header><section class="modal-card-body"><p><label class="label">Theme</label><div class="select"><select id="documenter-themepicker"><option value="documenter-light">documenter-light</option><option value="documenter-dark">documenter-dark</option></select></div></p><hr/><p>This document was generated with <a href="https://github.com/JuliaDocs/Documenter.jl">Documenter.jl</a> on <span class="colophon-date" title="Monday 16 May 2022 19:38">Monday 16 May 2022</span>. Using Julia version 1.7.0.</p></section><footer class="modal-card-foot"></footer></div></div></div></body></html>