Reference
Binary DD
ElementaryFluxModes.DDBinary
— MethodDDBinary(N, K) -> Any
Implement the Double Description method in binary form. The input variables are: -N
: the stoichiometric matrix with only forward reactions. If any fluxes are fixed then N has the form [N1 w] where N1 is the stoichiometric matrix for non-fixed fluxes, and w is the columns of fixed fluxes multiplied by their flux values. -K
: initial nullspace of N, in the form [I;K*] Output: -R
: binary elementary flux modes
ElementaryFluxModes.check_adjacency
— Methodcheck_adjacency(i, j, R) -> Bool
Check adjacency of columns i and j in r by making sure that there exists no other extreme ray in R whose zero set is a superset of the intersection of the zero sets of ray i and ray j.
ElementaryFluxModes.get_efms
— Methodget_efms(N::Matrix{Float64}; tol) -> Any
Calculate elementary flux modes of a stoichiometric matrix. Input: N
: the stoichiometric matrix of a network with only forward reactions. Output: matrix of the EFMs, each column is a different EFM, rows correspond to the reaction indices of the input N
.
ElementaryFluxModes.get_ofms
— Methodget_ofms(
N::Matrix{Float64},
fixed_fluxes::Vector{Int64},
flux_values::Vector{Float64}
) -> Any
Calculate the optimal flux modes of a pruned optimal solution. Arguments:
N
: the stoichiometric matrix of a pruned modelfixed_fluxes
: the indexes of the fixed fluxesflux_values
: the optimal values of the fixed fluxes
Output: matrix of the OFMs, each column is a different OFM, rows correspond to the reaction indices of the input N
.
ElementaryFluxModes.make_bitmap
— Methodmake_bitmap(row) -> Vector{Bool}
Return a bitmap of the input vector, where entries are 1 if the vector entry is greater than zero, and zero if the entries are zero. Return an error if any entries of the vector are less than zero.
ElementaryFluxModes.rational_nullspace
— Methodrational_nullspace(
A::Matrix;
tol
) -> Tuple{Matrix{Float64}, Any}
Helper function to calculate a nullspace of the matrix A, with all rational entries.
ElementaryFluxModes.zero_set
— Methodzero_set(vec) -> Vector
Return the zero set of a vector
Differentiating
ElementaryFluxModes._cost_matrix
— Method_cost_matrix(
EFMs::Vector{Dict{String, Float64}},
rid_pid,
rid_gcounts,
capacity::Vector{Tuple{String, Vector{String}, Float64}},
gene_product_molar_masses::Dict{String, Float64};
evaluate,
parameter_values
) -> Matrix{Any}
Calculate a matrix of the cost vectors of each EFM to each constraint. Entry (i,j) gives the total cost in constraint i to produce one unit objective flux through EFM j. Cost is calculated as ∑w(i)V(j)/kcat, where the variables are:
- 'w(i)': fraction of the ith enzyme pool that one mole of the enzyme uses up
- 'V(j)': flux through the reaction in EFM j
- 'kcat': turnover number of the enzyme.
ElementaryFluxModes.differentiate_efm
— Methoddifferentiate_efm(
EFMs::Vector{Dict{String, Float64}},
parameters::Vector{FastDifferentiation.Node},
rid_pid::Dict{String, FastDifferentiation.Node},
parameter_values::Dict{Symbol, Float64},
rid_gcounts::Dict{String, Dict{String, Float64}},
capacity::Vector{Tuple{String, Vector{String}, Float64}},
gene_product_molar_masses::Dict{String, Float64},
optimizer
) -> Any
Differentiate the usage of EFMs with respect to changes in the model parameters using implicit differentiation of the Lagrangian. Calculating the proportion of each EFM used in the optimal solution can be turned into the following LP: maximise ∑xᵢ subject to Dx = 1 where D is the cost matrix associated to each EFM in each enzyme pool. We then differentiate Lagrangian of this LP to calculate the differential of x by parameters.
Arguments:
- 'EFMs': a vector of dictionaries of EFMs
- 'parameters': vector of the model parameters
rid_pid
: dict of reaction ids => parameter idsparameter_values
: dict of parameter symbol => float valuerid_gcounts
: dict of reaction id => gene product countscapacity
: enzyme capacity limit of the enzyme constrained modelgene_product_molar_masses
: dict of gene product geneproductmolar_masses
Output: Matrix Aᵢⱼ of the the differential of EFM i with respect to parameter j: d(EFMᵢ)/d(pⱼ)
Utilities
ElementaryFluxModes.remove_linearly_dep_rows_qr
— Methodremove_linearly_dep_rows_qr(A) -> Any
Helper function to remove linearly dependent rows of the matrix A, taken from DifferentiableMetabolism.jl, using QR decomposition.
ElementaryFluxModes.reorder_ns
— Methodreorder_ns(A::Matrix) -> Tuple{Matrix, Vector{Int64}}
Helper function to reorder the rows of the nullspace so that it is in the form [I; K].