Markov Chain Utilities
CTMC and DTMC analysis tools.
The mc module provides general tools for analyzing Markov chains, including
steady-state and transient analysis for both continuous-time and discrete-time
Markov chains.
Key function categories:
Steady-state analysis:
ctmc_solve_reducible(),dtmc_solve_reducible(),ctmc_stochcomp()Transient analysis:
ctmc_transient(),ctmc_uniformization()Simulation:
ctmc_simulate(),ctmc_rand()Aggregation methods:
ctmc_courtois(),ctmc_takahashi(),ctmc_kms()State-space generation:
ctmc_ssg(),ctmc_ssg_reachability()Generator matrices:
ctmc_makeinfgen(),ctmc_multi()
Markov Chain Analysis (line_solver.api.mc)
Markov Chain analysis algorithms.
Native Python implementations for continuous-time and discrete-time Markov chain analysis.
- Key algorithms:
ctmc_solve: CTMC steady-state distribution ctmc_transient: CTMC transient analysis ctmc_uniformization: Uniformization method for transient analysis ctmc_stochcomp: Stochastic complementation dtmc_solve: DTMC steady-state distribution
- ctmc_solve(Q)[source]
Solve for steady-state probabilities of a CTMC.
Computes the stationary distribution π by solving πQ = 0 with normalization constraint Σπ = 1.
Handles reducible CTMCs by decomposing into strongly connected components and solving each separately.
- ctmc_makeinfgen(Q)[source]
Convert a matrix into a valid infinitesimal generator for a CTMC.
An infinitesimal generator has: - Row sums equal to zero - Non-positive diagonal elements - Non-negative off-diagonal elements
- ctmc_transient(Q, initial_dist, time_points, method='expm')[source]
Compute transient probabilities of a CTMC.
Calculates time-dependent state probabilities π(t) for specified time points using matrix exponential methods.
- Parameters:
- Returns:
Transient probabilities at each time point. Shape: (len(time_points), n) if multiple times, (n,) if single time
- Return type:
- ctmc_uniformization(Q, lambda_rate=None)[source]
Uniformize CTMC generator matrix.
Converts CTMC to an equivalent uniformized discrete-time chain for numerical analysis and simulation purposes.
The uniformized DTMC has transition matrix P = I + Q/λ where λ is the uniformization rate (max exit rate).
- ctmc_randomization(Q, initial_dist, time_points, precision=1e-10)[source]
Compute CTMC transient probabilities using randomization.
Uses Jensen’s randomization method (uniformization) to compute transient probabilities by converting the CTMC to a uniformized DTMC.
This method is numerically stable and avoids matrix exponentials.
- Parameters:
- Returns:
Transient probabilities at each time point
- Return type:
- ctmc_stochcomp(Q, keep_states, eliminate_states=None)[source]
Compute stochastic complement of CTMC.
Reduces the CTMC by eliminating specified states while preserving the steady-state distribution restricted to the kept states.
- Parameters:
- Returns:
‘S’: Stochastic complement (reduced generator)
’Q11’: Submatrix for kept states
’Q12’: Transitions from kept to eliminated
’Q21’: Transitions from eliminated to kept
’Q22’: Submatrix for eliminated states
’T’: Transient contribution matrix
- Return type:
dict containing
- ctmc_timereverse(Q, pi=None)[source]
Compute time-reversed CTMC generator.
The time-reversed generator Q* has elements: Q*_{ij} = π_j * Q_{ji} / π_i
- ctmc_simulate(Q, initial_state, max_time, max_events=10000, seed=None)[source]
Simulate CTMC sample path using Gillespie algorithm.
Generates a realization of the continuous-time Markov chain using the next-reaction method.
- Parameters:
- Returns:
‘states’: Array of visited states
’times’: Array of transition times
’sojourn_times’: Time spent in each state
- Return type:
dict with
- ctmc_isfeasible(Q, tolerance=1e-10)[source]
Check if matrix is a valid CTMC infinitesimal generator.
Validates: - Off-diagonal elements are non-negative - Row sums are zero - Diagonal elements are non-positive
- ctmc_ssg(sn, options=None)[source]
Generate complete CTMC state space for a queueing network.
Creates all possible network states including those not reachable from the initial state. For open classes, a cutoff parameter limits the maximum population to keep state space finite.
The state space is aggregated to show per-station-class job counts.
- Parameters:
- Returns:
state_space: Complete state space matrix (rows=states, cols=state components)
state_space_aggr: Aggregated state space (rows=states, cols=stations*classes)
state_space_hashed: Hashed state indices for lookup
node_state_space: Dictionary of per-node state spaces
sn: Updated network structure with space field populated
- Return type:
CtmcSsgResult containing
References
MATLAB: matlab/src/api/mc/ctmc_ssg.m
- ctmc_ssg_reachability(sn, options=None)[source]
Generate reachable CTMC state space for a queueing network.
Creates only the states reachable from the initial state through valid transitions. This is more efficient than ctmc_ssg for networks with constrained reachability.
- Parameters:
- Returns:
state_space: Reachable state space matrix
state_space_aggr: Aggregated state space (per station-class)
state_space_hashed: Hashed state indices
node_state_space: Dictionary of per-node state spaces
sn: Updated network structure
- Return type:
CtmcSsgResult containing
References
MATLAB: matlab/src/api/mc/ctmc_ssg_reachability.m
- class CtmcSsgResult(state_space, state_space_aggr, state_space_hashed, node_state_space, sn)[source]
Bases:
objectResult from CTMC state space generation.
- dtmc_solve(P)[source]
Solve for steady-state probabilities of a DTMC.
Computes the stationary distribution π by solving π(P - I) = 0 with normalization constraint Σπ = 1.
This leverages the CTMC solver by treating (P - I) as an infinitesimal generator.
- dtmc_solve_reducible(P, pin=None)[source]
Solve reducible DTMCs with transient states.
Handles DTMCs with multiple recurrent classes and transient states by: 1. Decomposing into strongly connected components (SCCs) 2. Identifying recurrent vs transient SCCs 3. Computing limiting distribution considering absorption from transient states
For a reducible DTMC with a single transient SCC, this computes the limiting distribution when starting from the transient states (e.g., class switching networks where jobs start in a transient class).
- dtmc_makestochastic(A)[source]
Convert matrix to row-stochastic transition matrix.
Normalizes each row to sum to 1. Rows with zero sum are replaced with uniform distribution.
- dtmc_isfeasible(P, tolerance=1e-10)[source]
Check if matrix is a valid DTMC transition matrix.
Validates: - All elements are non-negative - All row sums equal 1
- dtmc_simulate(P, initial_state, num_steps, seed=None)[source]
Simulate DTMC sample path.
Generates a realization of the discrete-time Markov chain for a specified number of steps.
- dtmc_timereverse(P, pi=None)[source]
Compute time-reversed DTMC transition matrix.
The time-reversed chain has transition probabilities: P*_{ij} = π_j * P_{ji} / π_i
- dtmc_stochcomp(P, keep_states, eliminate_states=None)[source]
Compute stochastic complement of DTMC.
Reduces the DTMC by eliminating specified states while preserving the steady-state distribution restricted to the kept states.
- dtmc_transient(P, initial_dist, steps)[source]
Compute transient probabilities of a DTMC.
Calculates π(n) = π(0) * P^n for each step from 0 to steps.
- dtmc_hitting_time(P, target_states)[source]
Compute mean hitting times to target states.
Calculates the expected number of steps to reach any target state from each starting state.
- class CourtoisResult(p, Qperm, Qdec, eps, epsMAX, P, B, q)[source]
Bases:
objectResult of Courtois decomposition.
- class KMSResult(p, p_1, Qperm, eps, epsMAX, pcourt)[source]
Bases:
objectResult of KMS aggregation-disaggregation.
- class TakahashiResult(p, p_1, pcourt, Qperm, eps, epsMAX)[source]
Bases:
objectResult of Takahashi aggregation-disaggregation.
- ctmc_courtois(Q, MS, q=None)[source]
Courtois decomposition for near-completely decomposable CTMCs.
Decomposes a large CTMC into macrostates and computes approximate steady-state probabilities using hierarchical aggregation.
- Parameters:
- Returns:
CourtoisResult with approximate solution and diagnostics
- Return type:
References
Original MATLAB: matlab/src/api/mc/ctmc_courtois.m Courtois, “Decomposability: Queueing and Computer System Applications”, 1977
- ctmc_kms(Q, MS, numSteps=10)[source]
Koury-McAllister-Stewart aggregation-disaggregation method.
Iteratively refines the Courtois decomposition solution using aggregation and disaggregation steps.
- Parameters:
- Returns:
KMSResult with refined solution
- Return type:
References
Original MATLAB: matlab/src/api/mc/ctmc_kms.m Koury, McAllister, Stewart, “Iterative Methods for Computing Stationary Distributions of Nearly Completely Decomposable Markov Chains”, 1984
- ctmc_takahashi(Q, MS, numSteps=10)[source]
Takahashi’s aggregation-disaggregation method.
Iteratively refines the Courtois decomposition solution using a different aggregation-disaggregation scheme.
- Parameters:
- Returns:
TakahashiResult with refined solution
- Return type:
References
Original MATLAB: matlab/src/api/mc/ctmc_takahashi.m Takahashi, “A Lumping Method for Numerical Calculations of Stationary Distributions of Markov Chains”, 1975
Continuous-Time Markov Chains (line_solver.api.ctmc)
The ctmc module contains algorithms specifically for continuous-time Markov
chains (CTMCs), including steady-state solvers and transient analysis.
Discrete-Time Markov Chains (line_solver.api.dtmc)
The dtmc module provides algorithms for discrete-time Markov chains (DTMCs),
including steady-state and transient analysis.