CHMC

Public functions

MarkovWeightedEFMs.CHMC.Standard.steady_state_efm_distributionFunction
steady_state_efm_distribution(#
    S::Matrix{<:Integer},
    v::Vector{<:Real},
    I::Int64=1;
    solver::Symbol=:Direct,
    issparse::Bool=false
)

Enumerate the EFMs from stoichiometry matrix S and compute the steady state probabilities of each EFM according to the discrete-time, cycle-history Markov chain.

S is a fully-connected, unimolecular, m by n stoichiometry matrix with m metabolites and n reactions.

v is the n-length steady state flux vector associated with S.

I is the initial starting state for rooting the cycle-history Markov chain. The choice of initial starting state does not affect the steady state EFM probabilities. The default is 1 and must be a whole number between 1:m.

solver is the type used for eigenvector calculations. Default is :Direct but may use iterative methods :IterativeSolver_gmres or :Arnoldi

issparse is true will use a sparse transition probability for the CHMC to save memory.

Example

julia> S = [#
 -1  0  0  0  0  0  0  0  0  0  1
  1 -1  1 -1  0  0  0  0  0  0  0
  0  1 -1  0 -1  1  0  0  0  0  0
  0  0  0  1  0  0 -1  0  0  0  0
  0  0  0  0  1 -1  1 -1  1 -1  0
  0  0  0  0  0  0  0  0  0  1 -1
  0  0  0  0  0  0  0  1 -1  0  0
];
julia> v = [3, 2, 1, 2, 3, 2, 2, 1, 1, 3, 3];
julia> res = steady_state_efm_distribution(S, v);
julia> res.e # EFM state sequences
6-element Vector{Vector{Int64}}:
 [3, 2, 3]
 [3, 2, 4, 5, 3]
 [3, 5, 3]
 [6, 1, 2, 4, 5, 6]
 [7, 5, 7]
 [6, 1, 2, 3, 5, 6]

julia> res.p # EFM probabilities
6-element Vector{Float64}:
 0.10638297872340426
 0.0425531914893617
 0.25531914893617025
 0.1914893617021277
 0.14893617021276595
 0.25531914893617025

julia> res.w # EFM weights
6-element Vector{Float64}:
 0.7142857142857142
 0.2857142857142857
 1.7142857142857144
 1.2857142857142858
 0.9999999999999999
 1.7142857142857144
source
MarkovWeightedEFMs.CHMC.Standard.steady_state_efm_distributionFunction
steady_state_efm_distribution(#
    T::Matrix{<:Real},
    I::Int64=1;
    solver::Symbol=:Direct
)

Enumerate the EFMs from (right) transition probability matrix T whose rows sum to one, and compute the steady state probabilities of each EFM according to the discrete-time, cycle-history Markov chain.

T is the discrete-time transition probability matrix with probabilities proportional to the outgoing fluxes.

I is the initial starting state for rooting the cycle-history Markov chain. The choice of initial starting state does not affect the steady state EFM probabilities. The default is 1 and must be a whole number between 1:m.

solver is the type used for eigenvector calculations. Default is :Direct but may use iterative methods :IterativeSolver_gmres or :Arnoldi.

Example

julia> T = [#
  0.0  1.0   0.0       0.0  0.0   0.0  0.0
  0.0  0.0   0.5       0.5  0.0   0.0  0.0
  0.0  0.25  0.0       0.0  0.75  0.0  0.0
  0.0  0.0   0.0       0.0  1.0   0.0  0.0
  0.0  0.0   0.333333  0.0  0.0   0.5  0.166667
  1.0  0.0   0.0       0.0  0.0   0.0  0.0
  0.0  0.0   0.0       0.0  1.0   0.0  0.0
];
julia> res = steady_state_efm_distribution(T);
julia> res.e # EFM state sequences
6-element Vector{Vector{Int64}}:
 [3, 2, 3]
 [3, 2, 4, 5, 3]
 [3, 5, 3]
 [6, 1, 2, 3, 5, 6]
 [7, 5, 7]
 [6, 1, 2, 4, 5, 6]

julia> res.p # EFM probabilities
6-element Vector{Float64}:
 0.13723110896294213
 0.035797649943428565
 0.26989203316869914
 0.20302350628312624
 0.13926980198123254
 0.2147858996605714
source
MarkovWeightedEFMs.CHMC.Standard.stoichiometry_to_transition_matrixFunction
stoichiometry_to_transition_matrix(S::Matrix{<:Real}, v::Vector{<:Real})

Convert stoichiometry matrix S with vector of steady state fluxes to a right stochastic transition probability matrix with rows summing to one.

S is the m by n stoichiometry matrix with m metabolites and n reactions.

v is the steady state flux vector with length n.

Examples

julia> S = [#
  -1  0  0  0  0  0  0  0  0  0  1
   1 -1  1 -1  0  0  0  0  0  0  0
   0  1 -1  0 -1  1  0  0  0  0  0
   0  0  0  1  0  0 -1  0  0  0  0
   0  0  0  0  1 -1  1 -1 -1  1  0
   0  0  0  0  0  0  0  1  0  0 -1
   0  0  0  0  0  0  0  0  1 -1  0
]
julia> v = [2, 2, 2, 2, 2, 2, 2, 2, 4]
julia> stoich_to_transition(S, v)
7x7 Matrix{Float64}:
 0.0  1.0   0.0       0.0  0.0   0.0  0.0
 0.0  0.0   0.5       0.5  0.0   0.0  0.0
 0.0  0.25  0.0       0.0  0.75  0.0  0.0
 0.0  0.0   0.0       0.0  1.0   0.0  0.0
 0.0  0.0   0.333333  0.0  0.0   0.5 0.166667
 1.0  0.0   0.0       0.0  0.0   0.0  0.0
 0.0  0.0   0.0       0.0  1.0   0.0  0.0
source
MarkovWeightedEFMs.CHMC.Standard.reshape_efm_matrixFunction
reshape_efm_matrix(ϕ::Matrix{Int64}, S::Matrix{<:Real})

Convert a matrix of EFMs ϕ to a nested vector of EFMs from a stoichiometry matrix S. Stoichiometry matrix may only contain unimolecular reactions.

ϕ is the n by k EFM matrix with n reactions (rows) and k EFMs (cols).

S is the m by n stoichiometry matrix with m metabolites (rows) and n reactions (cols).

Examples

julia> ϕ = [#
  1  1  0  0  0  0
  1  0  1  0  0  0
  0  0  1  0  0  1
  0  1  0  0  0  1
  1  0  0  1  0  0
  0  0  0  1  0  1
  0  1  0  0  0  1
  1  1  0  0  0  0
  0  0  0  0  1  0
  0  0  0  0  1  0
  1  1  0  0  0  0
]
julia> S = [#
  -1  0  0  0  0  0  0  0  0  0  1
   1 -1  1 -1  0  0  0  0  0  0  0
   0  1 -1  0 -1  1  0  0  0  0  0
   0  0  0  1  0  0 -1  0  0  0  0
   0  0  0  0  1 -1  1 -1 -1  1  0
   0  0  0  0  0  0  0  1  0  0 -1
   0  0  0  0  0  0  0  0  1 -1  0
]
julia> efm_vector = reshape_efm_matrix(ϕ, S)
6-element Vector{Vector{Int64}}:
 [1, 2, 3, 5, 6, 1]
 [1, 2, 4, 5, 6, 1]
 [2, 3, 2]
 [3, 5, 3]
 [5, 7, 5]
 [3, 2, 4, 5, 3]
source
MarkovWeightedEFMs.CHMC.Standard.reshape_efm_vectorFunction
reshape_efm_vector(ϕ::Vector{Vector{Int64}}, S::Matrix{<:Real})

Convert nested vector of EFM indices ϕ with length k to an n by k matrix of EFMs based on m by n stoichiometry matrix S. Stoichiometry matrix may only contain unimolecular reactions.

ϕ is the nested vector of EFMs with length k and elements corresponding to EFM metabolite indices in S.

S is the m by n stoichiometry matrix with m metabolites (rows) and n reactions (cols).

Examples

julia> ϕ = [#
  [1, 2, 3, 5, 6, 1],
  [1, 2, 4, 5, 6, 1],
  [2, 3, 2], [3, 5, 3], [5, 7, 5], [2, 4, 5, 3, 2]
]
julia> S = [#
  -1  0  0  0  0  0  0  0  0  0  1
   1 -1  1 -1  0  0  0  0  0  0  0
   0  1 -1  0 -1  1  0  0  0  0  0
   0  0  0  1  0  0 -1  0  0  0  0
   0  0  0  0  1 -1  1 -1 -1  1  0
   0  0  0  0  0  0  0  1  0  0 -1
   0  0  0  0  0  0  0  0  1 -1  0
]
julia> efm_matrix = reshape_efm_vector(ϕ, S)
11x6 Matrix{Int64}:
 1  1  0  0  0  0
 1  0  1  0  0  0
 0  0  1  0  0  1
 0  1  0  0  0  1
 1  0  0  1  0  0
 0  0  0  1  0  1
 0  1  0  0  0  1
 1  1  0  0  0  0
 0  0  0  0  1  0
 0  0  0  0  1  0
 1  1  0  0  0  0
source

Index