Front-end user interface
Flux balance analysis
COBREXA.flux_balance_analysis
— Methodflux_balance_analysis(
model::AbstractFBCModels.AbstractFBCModel;
kwargs...
)
Compute an optimal objective-optimizing solution of the given model
.
Most arguments are forwarded to optimized_values
.
Returns a tree with the optimization solution of the same shape as given by flux_balance_constraints
.
COBREXA.flux_balance_constraints
— Methodflux_balance_constraints(
model::AbstractFBCModels.AbstractFBCModel;
interface,
interface_name
) -> Any
A constraint tree that models the content of the given instance of AbstractFBCModel
.
The constructed tree contains subtrees fluxes
(with the reaction-defining "variables") and flux_stoichiometry
(with the metabolite-balance-defining constraints), and a single constraint objective
thad describes the objective function of the model.
Optionally if interface
is specified, an "interface block" will be created within the constraint tree for later use as a "module" in creating bigger models (such as communities) using interface_constraints
. The possible parameter values include:
nothing
– default, no interface is created:sbo
– the interface gets created from model's SBO annotations):identifier_prefixes
– the interface is guesstimated from commonly occurring adhoc reaction ID prefixes used in contemporary models:boundary
– the interface is created from all reactions that either only consume or only produce metabolites
Output interface name can be set via interface_name
.
See Configuration
for fine-tuning the default interface creation.
Flux variability analysis
COBREXA.flux_variability_analysis
— Methodflux_variability_analysis(
model::AbstractFBCModels.AbstractFBCModel;
objective_bound,
reactions,
optimizer,
settings,
workers
)
Perform a Flux Variability Analysis (FVA) on the model
, and return a dictionary of flux ranges where the model is able to perform optimally.
The constraint system is constructed using flux_balance_constraints
, and the variability is examined on all reaction's fluxes, or on the subset given optionally in reaction_subset
(e.g., reaction_subset = ["PFK", "ACALD"]
). The optimality tolerance can be specified with objective_bound
using e.g. relative_tolerance_bound
or absolute_tolerance_bound
; the default is 99% relative tolerance.
Parameter workers
may be used to enable parallel or distributed processing; the execution defaults to all available workers. Other parameters (esp. optimizer
) are internally forwarded to optimized_values
.
Use constraints_variability
to customize the FVA execution.
Parsimonious flux balance analysis
COBREXA.linear_parsimonious_flux_balance_analysis
— Methodlinear_parsimonious_flux_balance_analysis(
model::AbstractFBCModels.AbstractFBCModel;
tolerances,
kwargs...
)
Like parsimonious_flux_balance_analysis
, but uses the L1-metric parsimonious system given by linear_parsimonious_flux_balance_constraints
.
In turn, the solution is often faster, does not require a solver capable of quadratic objectives, and has many beneficial properties of the usual parsimonious solutions (such as the general lack of unnecessary loops). On the other hand, like with plain flux balance analysis there is no strong guarantee of uniqueness of the solution.
Solver configuration arguments are forwarded to parsimonious_optimized_values
.
COBREXA.linear_parsimonious_flux_balance_constraints
— Methodlinear_parsimonious_flux_balance_constraints(
model::AbstractFBCModels.AbstractFBCModel;
kwargs...
) -> ConstraintTrees.Tree{ConstraintTrees.Constraint}
Like parsimonious_flux_balance_constraints
, but uses a L1 metric for solving the parsimonious problem. The parsimonious_objective
constraint is thus linear.
Keyword arguments are forwarded to flux_balance_constraints
.
COBREXA.parsimonious_flux_balance_analysis
— Methodparsimonious_flux_balance_analysis(
model::AbstractFBCModels.AbstractFBCModel;
tolerances,
kwargs...
)
Compute a parsimonious flux solution for the model
, using the constraints given by parsimonious_flux_balance_constraints
.
In short, the objective value of the parsimonious solution should be the same as the one from flux_balance_analysis
, except the squared sum of reaction fluxes is minimized. If there are multiple possible fluxes that achieve a given objective value, parsimonious flux thus represents the "minimum energy" one, which is arguably more realistic.
Solver configuration arguments are forwarded to parsimonious_optimized_values
.
COBREXA.parsimonious_flux_balance_constraints
— Methodparsimonious_flux_balance_constraints(
model::AbstractFBCModels.AbstractFBCModel;
kwargs...
) -> Any
A constraint system like from flux_balance_constraints
, but with the parsimonious objective present under key parsimonious_objective
. Best used via parsimonious_flux_balance_analysis
.
Keyword arguments are forwarded to flux_balance_constraints
.
Minimization of metabolic adjustment
COBREXA.linear_metabolic_adjustment_minimization_analysis
— Methodlinear_metabolic_adjustment_minimization_analysis(
model::AbstractFBCModels.AbstractFBCModel,
args...;
optimizer,
settings,
reference_parsimonious_optimizer,
reference_parsimonious_settings,
reference_optimizer,
reference_settings,
kwargs...
)
Perform a linear minimization of metabolic adjustment analysis (l-MOMA) on model
. The reference is given by the second argument, which is either a reference_flux
or a reference_model
(the second argument is forwarded to linear_metabolic_adjustment_minimization_constraints
).
While the solution is "less uniquely defined" than with fully quadratic metabolic_adjustment_minimization_analysis
, the linear variant typically produces a sufficiently good result with much less resources. See documentation of linear_parsimonious_flux_balance_analysis
for some of the considerations.
COBREXA.linear_metabolic_adjustment_minimization_constraints
— Methodlinear_metabolic_adjustment_minimization_constraints(
model::AbstractFBCModels.AbstractFBCModel,
reference_model::AbstractFBCModels.AbstractFBCModel;
kwargs...
) -> Union{Nothing, ConstraintTrees.Tree{ConstraintTrees.Constraint}}
Like metabolic_adjustment_minimization_constraints
but the output constraints optimize the L1 distance from the linear-parsimonious solution of the reference_model
.
COBREXA.linear_metabolic_adjustment_minimization_constraints
— Methodlinear_metabolic_adjustment_minimization_constraints(
model::AbstractFBCModels.AbstractFBCModel,
reference_fluxes::ConstraintTrees.Tree;
_...
) -> ConstraintTrees.Tree{ConstraintTrees.Constraint}
Like metabolic_adjustment_minimization_constraints
but optimizes the L1 distance from reference_fluxes
.
Keyword arguments are discarded for compatibility with the other overload.
COBREXA.metabolic_adjustment_minimization_analysis
— Methodmetabolic_adjustment_minimization_analysis(
model::AbstractFBCModels.AbstractFBCModel,
args...;
optimizer,
settings,
reference_parsimonious_optimizer,
reference_parsimonious_settings,
reference_optimizer,
reference_settings,
kwargs...
)
Find a solution of the "minimization of metabolic adjustment" (MOMA) analysis for the model
, which is the "closest" feasible solution to the solution given in the second argument, which is either reference_fluxes
or reference_model
(see documentation of metabolic_adjustment_minimization_constraints
), in the sense of squared-sum distance. The minimized squared distance (the objective) is present in the result tree as minimal_adjustment_objective
.
If the second argument is a reference model, it is solved using a parsimonious_flux_balance_analysis
with the optimizer and settings parameters for the 2 steps set by keyword arguments prefixed by reference_
.
This is often used for models with smaller feasible region than the reference models (typically handicapped by a knockout, a nutritional deficiency or a similar perturbation). MOMA solution then gives an expectable "easiest" adjustment of the organism towards a somewhat working state.
Reference fluxes that do not exist in the model
are ignored (internally, the objective is constructed via squared_sum_error_value
).
COBREXA.metabolic_adjustment_minimization_constraints
— Methodmetabolic_adjustment_minimization_constraints(
model::AbstractFBCModels.AbstractFBCModel,
reference_model::AbstractFBCModels.AbstractFBCModel;
kwargs...
) -> Any
A slightly easier-to-use version of metabolic_adjustment_minimization_constraints
that computes the reference flux as the parsimonious optimal solution of the reference_model
. The reference flux is calculated using reference_optimizer
and reference_modifications
, which default to the optimizer
and settings
.
Other arguments are forwarded to the internal call of parsimonious_optimized_values
.
Returns nothing
if no feasible solution is found.
COBREXA.metabolic_adjustment_minimization_constraints
— Methodmetabolic_adjustment_minimization_constraints(
model::AbstractFBCModels.AbstractFBCModel,
reference_fluxes::ConstraintTrees.Tree;
_...
) -> Any
Keyword arguments are discarded for compatibility with the other overload.
Constraint systems for metabolite concentrations
COBREXA.log_concentration_constraints
— Methodlog_concentration_constraints(
model::AbstractFBCModels.AbstractFBCModel;
reactions,
metabolites,
metabolite_concentration_bound,
reaction_concentration_bound
) -> ConstraintTrees.Tree{ConstraintTrees.Constraint}
Build log-concentration-stoichiometry constraints for the model
, as used e.g. by max_min_driving_force_analysis
.
The output constraint tree contains a log-concentration variable for each metabolite
from model
, in subtree log_concentrations
. The total reactant log-concentrations for each reaction
are constrained in subtree log_concentration_stoichiometry
. By default, all reactions and metabolites in model
are included.
A concentration bound is given by parameter function concentration_bound
for each metabolite ID (the string ID is the single argument of the function); by default the function returns nothing
and no bounds are installed. The same is used for reactions with reaction_concentration_bound
.
Enzyme-mass-constrained models
COBREXA.Isozyme
— Typemutable struct Isozyme
A simple struct storing information about the isozyme composition, including subunit stoichiometry and turnover numbers. Use with enzyme_constrained_flux_balance_analysis
.
Fields
gene_product_stoichiometry::Dict{String, Float64}
: Mapping of gene product identifiers ("genes" in FBC model nomenclature) to their relative amount required to construct one unit of the isozyme.
kcat_forward::Union{Nothing, Float64}
: Turnover number for this isozyme catalyzing the forward direction of the reaction.kcat_reverse::Union{Nothing, Float64}
: Turnover number for this isozyme catalyzing the reverse direction of the reaction.
COBREXA.enzyme_constrained_flux_balance_analysis
— Methodenzyme_constrained_flux_balance_analysis(
model::AbstractFBCModels.AbstractFBCModel;
kwargs...
)
Perform the enzyme-constrained flux balance analysis on the model
and return the solved constraint system.
Arguments are forwarded to enzyme_constrained_flux_balance_constraints
; solver configuration arguments are forwarded to optimized_values
.
COBREXA.enzyme_constrained_flux_balance_constraints
— Methodenzyme_constrained_flux_balance_constraints(
model::AbstractFBCModels.AbstractFBCModel;
reaction_isozymes,
gene_product_molar_masses,
capacity,
interface,
interface_name
)
Construct a enzyme-constrained flux-balance constraint system, following the method in GECKO algorithm (refer to: Sánchez, Benjamín J., et al. "Improving the phenotype predictions of a yeast genome‐scale metabolic model by incorporating enzymatic constraints." Molecular systems biology 13.8 (2017): 935).
The enzyme mass constraints depend primarily on the available isozymes, given in parameter reaction_isozymes
, which is a mapping of reaction identifiers to descriptions of Isozyme
s that may catalyze the particular reactions. The isozymes are built from gene products, the mass of which is specified by gene_product_molar_masses
. In total, the amount of gene product building material is limited by capacity
.
capacity
may be a single number, which sets the mass limit for "all described enzymes". Alternatively, capacity
may be a vector of identifier-genes-limit triples that together form a constraint (identified by the given identifier) that limits the total sum of the listed genes to the given limit.
interface
and interface_name
are forwarded to flux_balance_constraints
.
COBREXA.simplified_enzyme_constrained_flux_balance_analysis
— Methodsimplified_enzyme_constrained_flux_balance_analysis(
model::AbstractFBCModels.AbstractFBCModel;
kwargs...
)
Perform the enzyme-constrained flux balance analysis on the model
and return the solved constraint system.
Arguments are forwarded to simplified_enzyme_constrained_flux_balance_constraints
; solver configuration arguments are forwarded to optimized_values
.
COBREXA.simplified_enzyme_constrained_flux_balance_constraints
— Methodsimplified_enzyme_constrained_flux_balance_constraints(
model;
reaction_isozymes,
gene_product_molar_masses,
capacity,
interface,
interface_name
)
Like enzyme_constrained_flux_balance_constraints
, but automatically selects a single "fastest" isozyme for each reaction direction. In turn, the system requires much less variables in the constraint system description, and usually solves more efficiently, for the price of possibly finding suboptimal solutions. The method follows the SMOMENT algorithm described in Bekiaris, P.S., Klamt, S. Automatic construction of metabolic models with enzyme constraints. BMC Bioinformatics 21, 19 (2020). https://doi.org/10.1186/s12859-019-3329-9.
Arguments are as with enzyme_constrained_flux_balance_constraints
, with a major difference in capacity
handling: the identifier lists (2nd elements of the triples given in the list) are not identifiers of gene products, but identifiers of reactions.
Community models
COBREXA.community_flux_balance_analysis
— Methodcommunity_flux_balance_analysis(args...; kwargs...)
Run the Community Flux Balance Analysis on several models. All arguments are forwarded to community_flux_balance_constraints
which constructs the model; this function returns the solved model.
COBREXA.community_flux_balance_constraints
— Methodcommunity_flux_balance_constraints(
model_abundances,
community_exchange_bounds;
interface,
interface_exchanges,
interface_biomass,
default_community_exchange_bound
) -> Any
Construct an instance of a linear community Flux Balance Analysis (cFBA) model. The relative abundances of the organisms are known in advance; this function predicts the maximum achievable community growth rate.
model_abundances
is a dictionary-like object that maps model identifiers to tuples of models (usually subtypes of AbstractFBCModel
) and their abundances (such as: "bug1" => (bug1, 0.5), ...
). community_exchange_bounds
is a dictionary-like object that can be additionally used to restrict selected community exchange reactions (keys should be reaction IDs, the values are converted to ConstraintTrees
-like bounds). Bounds otherwise default to parameter default_community_exchange_bound
, which itself defaults to nothing
(i.e., unbounded).
If required, constraint trees may be supplied instead of AbstracFBCModel
s in model_abundances
. These must provide an interface compatible with interface_exchanges
and interface_biomass
.
interface
is forwarded to flux_balance_constraints
. interface_exchanges
and interface_biomass
are used to pick up the correct interface part to contribute to the community exchanges and community biomass.
Production envelopes
COBREXA.objective_production_envelope
— Methodobjective_production_envelope(
model::AbstractFBCModels.AbstractFBCModel,
reactions::Vector{String};
breaks,
optimizer,
settings,
workers
)
Find the objective production envelope of the model
in the dimensions given by reactions
.
This runs a variability analysis of the fluxes given by flux_balance_constraints
to determine an applicable range for the dimensions, then splits the dimensions to equal-sized breaks (of count breaks
for each dimension, i.e. total breaks ^ length(reactions)
individual "multidimensional breaks") thus forming a grid, and returns an array of fluxes through the model objective with the individual reactions fixed to flux as given by the grid.
optimizer
and settings
are used to construct the optimization models.
The computation is distributed to the specified workers
, defaulting to all available workers.
Gap-filling
COBREXA.gap_filling_analysis
— Methodgap_filling_analysis(args...; kwargs...)
Run the gap-filling analysis on a constraint system specified by gap_filling_constraints
.
COBREXA.gap_filling_constraints
— Methodgap_filling_constraints(
model::AbstractFBCModels.AbstractFBCModel,
universal_model::AbstractFBCModels.AbstractFBCModel,
objective_target::Float64;
universal_reaction_cost,
max_cost,
known_fills,
kwargs...
) -> ConstraintTrees.Tree{ConstraintTrees.Constraint}
Make a gap-filling system from a FBC model with gaps and an universal FBC model that contains reactions to be added into the original model.
The output system will be constrainted to reach at least objective_target
flux through the objective function. Generally, this should be set to an arbitrary small value such as 0.05
.
universal_reaction_cost
should assign a numeric cost of inclusion of each of the reactions in the universal_model
; by default all are assigned equal weight of 1
. max_cost
puts an optional maximum limit on the cost, which may help the solver to avoid exploring unnecessarily complex solutions. known_fills
may contain previous solutions of the same system; these will be made infeasible in the output constraint system in order to allow discovery of other ones.
Additional arguments are forwarded to flux_balance_constraints
that converts model
to constraints.
COBREXA.gap_filling_constraints
— Methodgap_filling_constraints(
;
system,
stoichiometry,
universal_fluxes,
universal_stoichiometry,
flux_cost,
max_cost,
known_fills
)
Make a gap-fillign system from a non-solving system
with a separated stoichiometry
, filling in possible fluxes from universal_fluxes
that are balanced with universal_stoichiometry
flux_cost
can be used to assign a weight to each given universal flux; the total cost is bounded by max_cost
.
known_fills
may contain previous solutions to be avoided; these are processed by gap_filling_known_fill_constraint
and attached to the system.
stoichiometry
needs to be extended to construct the final constraints, thus it should not be present in the original system
.
COBREXA.gap_filling_known_fill_constraint
— Methodgap_filling_known_fill_constraint(
fill_flags::ConstraintTrees.Tree{ConstraintTrees.Constraint},
known_flags::ConstraintTrees.Tree{Float64}
) -> ConstraintTrees.Constraint
Produce a constraint that can be added to the system made by gap_filling_constraints
to avoid repeating of a solution that includes reactions already generated by another solution, as represented by the solved fill_flags
.
Parameter fill_flags
are the gapfilling flags of the given constraint system, parameter known_flags
is expected to contain the solved fill_flags
for the solution that is to be avoided.
Knockout models
COBREXA.gene_knockout_constraints
— Methodgene_knockout_constraints(
fluxes::ConstraintTrees.Tree{ConstraintTrees.Constraint},
knockout_genes,
model::AbstractFBCModels.AbstractFBCModel
) -> ConstraintTrees.Tree{ConstraintTrees.Constraint}
Make a ConstraintTree
that simulates a gene knockout of knockout_genes
in the model
and disables corresponding fluxes
accordingly.
Keys of the fluxes must correspond to the reaction identifiers in the model
.
knockout_genes
may be any collection that support element tests using in
. Since the test is done many times, a Set
is a preferred contained for longer lists of genes.
All constraints are equality constraints returned in a single flat ConstraintTree
.
COBREXA.gene_knockout_constraints
— Methodgene_knockout_constraints(
fluxes::ConstraintTrees.Tree{ConstraintTrees.Constraint},
knockout_gene::String,
model::AbstractFBCModels.AbstractFBCModel
) -> ConstraintTrees.Tree{ConstraintTrees.Constraint}
Convenience overload of gene_knockout_constraints
for knocking out a single gene (without the necessity to store the gene identifier in a singleton container).
COBREXA.gene_knockouts
— Functiongene_knockouts(model::AbstractFBCModels.AbstractFBCModel)
gene_knockouts(
model::AbstractFBCModels.AbstractFBCModel,
gene_combinations::Vector{<:Union{String, Tuple{Vararg{String, N}} where N}};
kwargs...
)
Compute the objective value of the model
for all knockouts specified by gene_combinations
, which is a vector of gene IDs or tuples of gene IDs that are knocked out in groups.
Returns a vector in the same order as gene_combinations
.
Extra arguments (mainly, the optimizer
) are forwarded to screen_optimization_model
.
Loopless flux balance analysis
COBREXA.loopless_flux_balance_analysis
— Methodloopless_flux_balance_analysis(
model::AbstractFBCModels.AbstractFBCModel;
kwargs...
)
Perform the loopless flux balance analysis on the model
, returning the solved constraint system.
Arguments are forwarded to loopless_flux_balance_constraints
(see the documentation for the description of the constraint system); solver configuration arguments are forwarded to optimized_values
.
COBREXA.loopless_flux_balance_constraints
— Methodloopless_flux_balance_constraints(
model::AbstractFBCModels.AbstractFBCModel;
flux_infinity_bound,
driving_force_nonzero_bound,
driving_force_infinity_bound
) -> ConstraintTrees.Tree{ConstraintTrees.Constraint}
Construct a flux-balance constraint system from model
with added quasi-thermodynamic constraints that ensure that thermodynamically infeasible internal cycles can not occur. The method is closer described by: Schellenberger, Lewis, and, Palsson. "Elimination of thermodynamically infeasible loops in steady-state metabolic models.", Biophysical journal, 2011`.
The loopless condition comes with a performance overhead: the computation needs to find the null space of the stoichiometry matrix (essentially inverting it); and the subsequently created optimization problem contains binary variables for each internal reaction, thus requiring a MILP solver and a potentially exponential solving time.
Internally, the system is constructed by combining flux_balance_constraints
and loopless_constraints
.
The arguments driving_force_max_bound
and driving_force_nonzero_bound
set the bounds (possibly negated ones) on the virtual "driving forces" (G_i in the paper).
Max-Min Driving Force analysis
COBREXA.max_min_driving_force_analysis
— Methodmax_min_driving_force_analysis(
model::AbstractFBCModels.AbstractFBCModel;
kwargs...
)
Perform the max-min driving force analysis on the model
, returning the solved constraint system.
Arguments are forwarded to max_min_driving_force_constraints
(see the documentation for the description of the constraint system); solver configuration arguments are forwarded to optimized_values
.
COBREXA.max_min_driving_force_constraints
— Methodmax_min_driving_force_constraints(
model::AbstractFBCModels.AbstractFBCModel;
reaction_standard_gibbs_free_energies,
reference_flux,
concentration_ratios,
constant_concentrations,
ignored_metabolites,
proton_metabolites,
water_metabolites,
concentration_lower_bound,
concentration_upper_bound,
T,
R,
reference_flux_atol
)
Create max-min driving force constraint system from model
, using the supplied reaction standard Gibbs free energies in reaction_standard_gibbs_free_energies
.
The method is described by: Noor, et al., "Pathway thermodynamics highlights kinetic obstacles in central metabolism.", PLoS computational biology, 2014.
reference_flux
sets the directions of each reaction in model
. The scale of the values is not important, only the direction is examined (w.r.t. reference_flux_atol
tolerance). Ideally, the supplied reference_flux
should be completely free of internal cycles, which enables the thermodynamic consistency. To get the cycle-free flux, you can use loopless_flux_balance_analysis
(computationally demanding, but gives thermodynamically consistent solutions), parsimonious_flux_balance_analysis
or linear_parsimonious_flux_balance_analysis
(which is computationally simple, but the consistency is not guaranteed).
Internally, log_concentration_constraints
is used to lay out the base structure of the problem.
Following arguments are set optionally:
water_metabolites
,proton_metabolites
andignored_metabolites
allow to completely ignore constraints on a part of the metabolite set, which is explicitly recommended especially for water and protons (since the analyses generally assume aqueous environment of constant pH)constant_concentrations
can be used to fix the concentrations of the metabolitesconcentration_lower_bound
andconcentration_upper_bound
set the default concentration bounds for all other metabolitesconcentration ratios
is a dictionary that assigns a tuple of metabolite-metabolite-concentration ratio constraint names; e.g. ATP/ADP ratio can be fixed to five-times-more-ATP by settingconcentration_ratios = Dict("adenosin_ratio" => ("atp", "adp", 5.0))
T
andR
default to the usual required thermodynamic constraints in the expected units (the defaults assume the "usual" units, valuing 298.15 K and ~0.008314 kJ/K/mol, respectively). These multiply the log-concentrations to obtain the actual Gibbs energies, and thus driving forces.
Sampling
COBREXA.flux_sample
— Methodflux_sample(
model::AbstractFBCModels.AbstractFBCModel;
seed,
objective_bound,
method,
n_chains,
collect_iterations,
optimizer,
settings,
workers,
kwargs...
)
Run a sampling algorithm on the near-optimal feasible space of the model
(as specified by objective_bound
). By default, the sampling algorithm is ACHR (the method
parameter is defaulted to sample_chain_achr
). The sampling algorithm is ran for n_chains
and the iterations for collecting the sampled values are specified by collect_iterations
.
optimizer
is used to generate the warmup (with settings
) for the sampler using the usual unidimensional maximum-variability fluxes (as from flux_variability_analysis
). All computations are parallelized across workers
.
Extra arguments are forwarded to sample_constraints
. Eventually the arguments will reach the method
function, so extra arguments can be also used to customize the methods (e.g., by setting the epsilon
for the ACHR sampler).