DifferentiableMetabolism
Documentation for DifferentiableMetabolism.
DifferentiableMetabolism.DifferentiableModel
DifferentiableMetabolism.Enzyme
DifferentiableMetabolism.ReferenceFunctions
DifferentiableMetabolism.ReferenceFunctions
Base.show
DifferentiableMetabolism._add_gecko_enzyme_variable_as_function
DifferentiableMetabolism._build_gecko_equality_enzyme_constraints
DifferentiableMetabolism._build_smoment_kcat_coupling
DifferentiableMetabolism._dg
DifferentiableMetabolism._differentiable_michaelis_menten_gecko_opt_problem
DifferentiableMetabolism._differentiable_thermodynamic_smoment_opt_problem
DifferentiableMetabolism._differentiate_kkt!
DifferentiableMetabolism._get_exponent_or_cutoff
DifferentiableMetabolism._linear_solver
DifferentiableMetabolism._make_differentiable_model
DifferentiableMetabolism._remove_lin_dep_rows
DifferentiableMetabolism._saturation
DifferentiableMetabolism._scale_derivatives
DifferentiableMetabolism._solve_model!
DifferentiableMetabolism.check_scaling
DifferentiableMetabolism.check_scaling
DifferentiableMetabolism.differentiate
DifferentiableMetabolism.differentiate!
DifferentiableMetabolism.isozyme_to_enzyme
DifferentiableMetabolism.make_analytic_var_derivative
DifferentiableMetabolism.make_derivatives
DifferentiableMetabolism.make_derivatives_with_equality_scaling
DifferentiableMetabolism.make_param_deriv_with_scaling
DifferentiableMetabolism.make_symbolic_param_derivative
DifferentiableMetabolism.make_var_deriv_with_scaling
DifferentiableMetabolism.prune_model
DifferentiableMetabolism.scaling_factor
DifferentiableMetabolism.update_E!
DifferentiableMetabolism.update_M!
DifferentiableMetabolism.update_Q!
DifferentiableMetabolism.update_c!
DifferentiableMetabolism.update_d!
DifferentiableMetabolism.update_h!
DifferentiableMetabolism.with_parameters
DifferentiableMetabolism.with_parameters
DifferentiableMetabolism.DifferentiableModel
— Typemutable struct DifferentiableModel
A struct representing a generalized differentiable constraint based model. Format of the model is:
minimize ½ * x' * Q(θ) * x + c(θ)' * x
s.t. E(θ) * x = d(θ)
M(θ) * x ≤ h(θ)
Nominally, every component can be differentiated.
Q::Function
c::Function
E::Function
d::Function
M::Function
h::Function
θ::Vector{Float64}
analytic_var_derivs::Function
analytic_par_derivs::Function
var_ids::Vector{String}
param_ids::Vector{String}
Note, to ensure differentiability, preprocessing of the model used to derive the DifferentiableModel
is required. In short, only an active solution may be differentiated, this requires that:
- the base model does not possess any isozymes, each reaction may be catalyzed by one enzyme (complex) only,
- all the reactions should be unidirectinal,
- the kcats in
rid_enzyme
are for the appropriate direction used in the model, - all rids in
rid_enzyme
are used in the model.
DifferentiableMetabolism.Enzyme
— Typemutable struct Enzyme
A struct used to store information about an enzyme. A simplified version of COBREXA.Isozyme
, which stores only the turnover number (kcat) that should be used by the differentiable model, as well as information about the protein complex (stoichiometry and molar mass).
kcat::Float64
gene_product_count::Dict{String, Int64}
gene_product_mass::Dict{String, Float64}
DifferentiableMetabolism.ReferenceFunctions
— Typestruct ReferenceFunctions
A struct used to store information the original optimization problem.
Q::Function
c::Function
E::Function
d::Function
M::Function
h::Function
nx::Int64
neq::Int64
nineq::Int64
nθ::Int64
DifferentiableMetabolism.ReferenceFunctions
— MethodReferenceFunctions(
dm::DifferentiableModel
) -> ReferenceFunctions
A constructor for ReferenceFunctions
.
Base.show
— MethodPretty printing of a differentiable model.
DifferentiableMetabolism._add_gecko_enzyme_variable_as_function
— Method_add_gecko_enzyme_variable_as_function(
rid_enzyme,
original_rid,
E_components,
col_idx,
gene_ids
)
Helper function to add an column into the enzyme stoichiometric matrix parametrically.
DifferentiableMetabolism._build_gecko_equality_enzyme_constraints
— Method_build_gecko_equality_enzyme_constraints(
gm::COBREXA.GeckoModel,
rid_enzyme
) -> Tuple{@NamedTuple{row_idxs::Vector{Int64}, col_idxs::Vector{Int64}, coeff_tuple::Vector{Tuple{Float64, String}}}, Vector}
Helper function to build the equality enzyme constraints.
DifferentiableMetabolism._build_smoment_kcat_coupling
— Method_build_smoment_kcat_coupling(
smm::COBREXA.SMomentModel,
rid_enzyme
) -> Tuple{Vector{Int64}, Vector}
Helper function to build kcat coupling in parametric form for smoment problems.
DifferentiableMetabolism._dg
— Method_dg(
model,
rid_enzyme,
rid_dg0,
rid,
mangled_rid,
θ;
RT,
ignore_reaction_ids,
ignore_metabolite_ids
) -> Any
Helper function to assign the thermodynamic driving force to a reaction.
DifferentiableMetabolism._differentiable_michaelis_menten_gecko_opt_problem
— Method_differentiable_michaelis_menten_gecko_opt_problem(
gm::COBREXA.GeckoModel,
rid_enzyme::Dict{String, Enzyme},
rid_dg0::Dict{String, Float64},
rid_km::Dict{String, Dict{String, Float64}};
ϵ,
atol,
digits,
RT,
ignore_reaction_ids,
ignore_metabolite_ids
) -> Tuple{DifferentiableMetabolism.var"#64#71"{COBREXA.GeckoModel}, DifferentiableMetabolism.var"#E#69"{DifferentiableMetabolism.var"#Se#67"{Float64, Vector{Any}, Vector{Any}, COBREXA.GeckoModel, Dict{String, Enzyme}, Dict{String, Float64}, Dict{String, Dict{String, Float64}}, Vector{_A}, @NamedTuple{row_idxs::Vector{Int64}, col_idxs::Vector{Int64}, coeff_tuple::Vector{Tuple{Float64, String}}}, Int64, Int64}, _A1, Int64} where {_A<:Tuple{String, Any, Float64}, _A1}, DifferentiableMetabolism.var"#63#70"{_A, Int64} where _A, DifferentiableMetabolism.var"#65#72"{_A, Int64} where _A, DifferentiableMetabolism.var"#66#73"{Vector{Float64}, Vector{Float64}}, Any}
Return optimization problem where thermodynamic and saturation effects are incorporated into the gecko problem, but in differentiable format.
DifferentiableMetabolism._differentiable_thermodynamic_smoment_opt_problem
— Method_differentiable_thermodynamic_smoment_opt_problem(
smm::COBREXA.SMomentModel,
rid_enzyme::Dict{String, Enzyme},
rid_dg0::Dict{String, Float64};
ϵ,
atol,
digits,
RT,
ignore_reaction_ids,
ignore_metabolite_ids
) -> Tuple{DifferentiableMetabolism.var"#86#91"{COBREXA.SMomentModel}, DifferentiableMetabolism.var"#84#89", DifferentiableMetabolism.var"#85#90", DifferentiableMetabolism.var"#M#94"{_A, DifferentiableMetabolism.var"#kcat_thermo_coupling#92"{Float64, Vector{Any}, Vector{Any}, COBREXA.SMomentModel, Dict{String, Enzyme}, Dict{String, Float64}, Vector{_A1}, Vector{Int64}}} where {_A, _A1<:Tuple{String, Any, Float64}}, DifferentiableMetabolism.var"#88#95"{COBREXA.SMomentModel, _A, _B, Vector{Float64}, Vector{Float64}} where {_A, _B}, Vector}
Return structures that will allow the most basic form of smoment to be solved. No enzyme constraints allowed. Most effective enzyme is the only GRR. Assume unidirectional reactions.
DifferentiableMetabolism._differentiate_kkt!
— Method_differentiate_kkt!(
x,
ν,
λ,
A,
B,
dx,
diffmodel::DifferentiableModel;
use_analytic
)
Implicitly differentiate a convex quadratic or linear program using the KKT conditions. Suppose F(z(θ), θ) = 0
represents the optimality (KKT) conditions of some optimization problem specified by diffmodel
, where z
is a vector of the primal, x
, and dual, ν
and λ
, solutions. Then, mutate A = ∂₁F(z(θ), θ)
, B = ∂₂F(z(θ), θ)
, and the optimal solution.
This function is called by differentiate!
.
DifferentiableMetabolism._get_exponent_or_cutoff
— Method_get_exponent_or_cutoff(x; atol) -> Any
Return log₁₀(|x|)
or x
if |x| ≤ atol
.
DifferentiableMetabolism._linear_solver
— Method_linear_solver(A, B, dx)
Separate linear solve to have the ability to swap it our for something faster.
DifferentiableMetabolism._make_differentiable_model
— Method_make_differentiable_model(
c,
_E,
_d,
_M,
_h,
θ,
var_ids,
param_ids;
scale_equality,
scale_inequality
) -> DifferentiableModel
Internal helper function to construct a basic linear program representing a differentiable metabolic model.
DifferentiableMetabolism._remove_lin_dep_rows
— Method_remove_lin_dep_rows(A; ϵ, atol, digits) -> Any
Remove linearly dependent rows in A
deduced by reducing the matrix to row echelon form. Adjust sensitivity to numerical issues with ϵ
. After row echelon reduction, all elements are rounded to digits
after the decimal. Rows of all zeros (absolute value of the each element in row ≤ atol
) are removed. Beware, this operation is expensive for very large matrices.
DifferentiableMetabolism._saturation
— Method_saturation(
model,
rid_enzyme,
rid_km,
rid,
mangled_rid,
θ;
ignore_reaction_ids,
ignore_metabolite_ids
) -> Any
A helper function to incorporate saturation effects.
DifferentiableMetabolism._scale_derivatives
— Method_scale_derivatives(x, dx, diffmodel::DifferentiableModel)
Scale the derivatives: dx/dy => dlog(x)/dlog(y)
. Note, only scales the primal variables, not the duals.
DifferentiableMetabolism._solve_model!
— Method_solve_model!(
x,
ν,
λ,
diffmodel::DifferentiableModel,
optimizer;
modifications
)
Solve the optimization problem, aka the forward pass.
This function is called by differentiate!
.
DifferentiableMetabolism.check_scaling
— Methodcheck_scaling(
mat;
atol
) -> Union{Nothing, Tuple{Union{Nothing, Tuple}, Union{Nothing, Tuple}, Union{Nothing, Tuple}}}
Return the best case (if rescaled using scaling_factor
) and current scaling of a matrix mat
. Scaling is defined as the largest difference between the exponent of the smallest and largest value in each row of a matrix. To set the cut-off for when a coefficient is considered to be zero, change atol
.
DifferentiableMetabolism.check_scaling
— Methodcheck_scaling(
diffmodel::DifferentiableModel;
atol,
verbose
) -> Union{Nothing, Tuple{Union{Nothing, Tuple}, Union{Nothing, Tuple}, Union{Nothing, Tuple}}}
Helper function to check the scaling of the equality constraint matrix.
DifferentiableMetabolism.differentiate!
— Methoddifferentiate!(
x,
ν,
λ,
A,
B,
dx,
diffmodel::DifferentiableModel,
optimizer;
modifications,
use_analytic,
scale_output
)
Differentiates diffmodel
using in-place functions as much as possible. The meaning of the variables are:
x
: primal variablesν
: equality dual variablesλ
: inequality dual variablesA
: variable derivativesB
: parameter derivitivesdx
: the sensitivities of the variablesdiffmodel
: the model to differantiateoptimizer
: the solver used in the forward passmodifications
: COBREXA functions forwarded to the solveruse_analytic
: use the analytic derivatives (need to be supplied)scale_output
: flag if the primal variable sensitivities should be scaled bydx/dy => dlog(x)/dlog(y)
DifferentiableMetabolism.differentiate
— Methoddifferentiate(
diffmodel::DifferentiableModel,
optimizer;
use_analytic,
scale_output,
modifications
) -> Tuple{Vector{Float64}, Union{Vector, Matrix{Float64}}}
Solve and differentiate an optimization problem using the optimality conditions. The output can be scaled relative to the parameters and the solved variables with scale_output
. Optimizer modifications (from COBREXA.jl) can be supplied through modifications
. Analytic derivatives of the optimality conditions can be used by setting use_analytic
to true.
Internally calls differentiate!
, the in-place variant of this function, constructs the arguments internally.
DifferentiableMetabolism.isozyme_to_enzyme
— Methodisozyme_to_enzyme(
isozyme::COBREXA.Isozyme,
gene_product_mass_lookup::Union{Dict{String, Float64}, Function};
direction
) -> Enzyme
Convert a COBREXA.Isozyme
to an Enzyme
using the turnover number in direction
(either :forward
or :reverse
), as well as gene_product_mass_lookup
, which is either a function or a dictionary mapping a gene product id to its molar mass.
DifferentiableMetabolism.make_analytic_var_derivative
— Methodmake_analytic_var_derivative(
dm::DifferentiableModel
) -> DifferentiableMetabolism.var"#103#104"{DifferentiableModel}
Creates the analytic derivatives of the variables.
DifferentiableMetabolism.make_derivatives
— Methodmake_derivatives(diffmodel::DifferentiableModel)
Creates analytic derivative functions of the diffmodel
internals, using symbolic variables. Note, this function can take some time to construct the derivatives, but substantially speeds up repeated calls to differentiate
.
See also make_derivatives_with_equality_scaling
for a variant that allows the equality constraints to be rescaled.
DifferentiableMetabolism.make_derivatives_with_equality_scaling
— Methodmake_derivatives_with_equality_scaling(
rf::ReferenceFunctions
) -> Tuple{DifferentiableMetabolism.var"#108#109"{ReferenceFunctions}, DifferentiableMetabolism.var"#105#107"{_A, SparseArrays.SparseMatrixCSC{Symbolics.Num, Int64}} where _A}
Returns analytic derivatives that allow the equality constraint to be scaled.
See also make_derivatives
.
DifferentiableMetabolism.make_param_deriv_with_scaling
— Methodmake_param_deriv_with_scaling(
rf::ReferenceFunctions
) -> DifferentiableMetabolism.var"#105#107"{_A, SparseArrays.SparseMatrixCSC{Symbolics.Num, Int64}} where _A
Return the parameter derivative that allows scaling of the equality constraint.
DifferentiableMetabolism.make_symbolic_param_derivative
— Methodmake_symbolic_param_derivative(
dm::DifferentiableModel
) -> DifferentiableMetabolism.var"#100#102"
Creates the analytic derivatives of the parameters using symbolic programming.
DifferentiableMetabolism.make_var_deriv_with_scaling
— Methodmake_var_deriv_with_scaling(
rf::ReferenceFunctions
) -> DifferentiableMetabolism.var"#108#109"{ReferenceFunctions}
Return the variable derivative that allows scaling of the equality constraint.
DifferentiableMetabolism.prune_model
— Methodprune_model(
model::COBREXA.StandardModel,
reaction_fluxes;
atol,
verbose
) -> COBREXA.StandardModel
Return a simplified version of model
that contains only reactions (and the associated genes and metabolites) that are active, i.e. carry fluxes (from reaction_fluxes
) absolutely bigger than atol
. All reactions are set unidirectional based on reaction_fluxes
.
DifferentiableMetabolism.scaling_factor
— Methodscaling_factor(A, b; atol) -> Any
Return scaling factor used to rescale constraints in matrix format (row-wise). Arguments are assumed to be of the form A ∝ b
where ∝
are usually =
or ≤
in the context of an optimization problem. To set the cut-off for when a coefficient is considered to be zero, change atol
.
See also check_scaling
.
DifferentiableMetabolism.update_E!
— Methodupdate_E!(diffmodel::DifferentiableModel, E) -> Any
Update E.
DifferentiableMetabolism.update_M!
— Methodupdate_M!(diffmodel::DifferentiableModel, M) -> Any
Update M.
DifferentiableMetabolism.update_Q!
— Methodupdate_Q!(diffmodel::DifferentiableModel, Q) -> Any
Update Q.
DifferentiableMetabolism.update_c!
— Methodupdate_c!(diffmodel::DifferentiableModel, c) -> Any
Update c.
DifferentiableMetabolism.update_d!
— Methodupdate_d!(diffmodel::DifferentiableModel, d) -> Any
Update d.
DifferentiableMetabolism.update_h!
— Methodupdate_h!(diffmodel::DifferentiableModel, h) -> Any
Update h.
DifferentiableMetabolism.with_parameters
— Methodwith_parameters(
gm::COBREXA.GeckoModel,
rid_enzyme::Dict{String, Enzyme};
rid_dg0,
rid_km,
mid_concentration,
scale_equality,
scale_inequality,
ϵ,
atol,
digits,
RT,
ignore_reaction_ids,
ignore_metabolite_ids
) -> DifferentiableModel
Construct a DifferentiableModel
from a COBREXA.GeckoModel
. Each variable in gm
is differentiated with respect to the kcats in the dictionary rid_enzyme
, which is a dictionary mapping reaction ids to Enzyme
s. Enzyme constraints are only taken with respect to the entries of rid_enzyme
. Optionally, incorporate thermodynamic and saturation constraints through rid_dg0
and rid_km
. If either thermodynamic or saturation (or both) constraints are added, then metabolite concentrations, mid_concentration
, need to be supplied as well.
Note, thermodynamic parameters require special attention. To ignore some reactions when calculating the thermodynamic factor, include them in ignore_reaction_ids
. The units used for the Gibbs free energy are kJ/mol, and concentrations should be in molar. Importantly, the standard Gibbs free energy of reaction is assumed to be given in the forward direction of the reaction in the model.
Internally, ϵ
, atol
, and digits
, are forwarded to _remove_lin_dep_rows
. Optionally, scale the equality constraints with scale_equality
, see scaling_factor
for more information.
DifferentiableMetabolism.with_parameters
— Methodwith_parameters(
smm::COBREXA.SMomentModel,
rid_enzyme::Dict{String, Enzyme};
rid_dg0,
mid_concentration,
scale_equality,
scale_inequality,
ϵ,
atol,
digits,
RT,
ignore_reaction_ids,
ignore_metabolite_ids
) -> DifferentiableModel
Construct a DifferentiableModel
from a COBREXA.SMomentModel
. Each variable in smm
is differentiated with respect to the kcats in the dictionary rid_enzyme
, which is a dictionary mapping reaction ids to Enzyme
s. Enzyme constraints are only taken with respect to the entries of rid_enzyme
. Optionally, include thermodynamic parameters. Add the standard Gibbs free energy of reactions with rid_dg0
, and the metabolite concentrations that are to be used as parameters with mid_concentration
(both arguments are dictionaries mapping reaction or metabolite ids to values).
Note, thermodynamic parameters require special attention. To ignore some reactions when calculating the thermodynamic factor, include them in ignore_reaction_ids
. The units used for the Gibbs free energy are kJ/mol, and concentrations should be in molar. Importantly, the standard Gibbs free energy of reaction is assumed to be given in the forward direction of the reaction in the model.
The analytic derivative of the optimality conditions with respect to the parameters can be supplied through analytic_parameter_derivatives
. Internally, ϵ
, atol
, and digits
, are forwarded to _remove_lin_dep_rows
.
Note, to ensure differentiability, preprocessing of the model is required. In short, only an active solution may be differentiated, this required that:
- the model does not possess any isozymes
- all the reactions should be unidirectinal
- the kcats in
rid_enzyme
are for the appropriate direction used in the model - all rids in
rid_enzyme
are used in the model