DifferentiableMetabolism

Documentation for DifferentiableMetabolism.

DifferentiableMetabolism.DifferentiableModelType
mutable 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.EnzymeType
mutable 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.ReferenceFunctionsType
struct 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

Base.showMethod

Pretty printing of a differentiable model.

DifferentiableMetabolism._build_gecko_equality_enzyme_constraintsMethod
_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._dgMethod
_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_problemMethod
_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_problemMethod
_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._make_differentiable_modelMethod
_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_rowsMethod
_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._saturationMethod
_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_derivativesMethod
_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.check_scalingMethod
check_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_scalingMethod
check_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!Method
differentiate!(
    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 variables
  • A: variable derivatives
  • B: parameter derivitives
  • dx: the sensitivities of the variables
  • diffmodel: the model to differantiate
  • optimizer: the solver used in the forward pass
  • modifications: COBREXA functions forwarded to the solver
  • use_analytic: use the analytic derivatives (need to be supplied)
  • scale_output: flag if the primal variable sensitivities should be scaled by dx/dy => dlog(x)/dlog(y)
DifferentiableMetabolism.differentiateMethod
differentiate(
    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_enzymeMethod
isozyme_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_derivatives_with_equality_scalingMethod
make_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_scalingMethod
make_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.prune_modelMethod
prune_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_factorMethod
scaling_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.with_parametersMethod
with_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 Enzymes. 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_parametersMethod
with_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 Enzymes. 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