LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
ctmc_solve_reducible_blkdecomp.m
1function [pi,pis,pi0,scc,isrec] = ctmc_solve_reducible_blkdecomp(Q,pi0,options)
2% [PI,PIS,PI0,SCC,ISREC] = CTMC_SOLVE_REDUCIBLE_BLKDECOMP(Q, PI0, OPTIONS)
3%
4% Compute limiting distribution for a CTMC with reducible generator Q
5% using direct block decomposition on the generator matrix.
6%
7% Algorithm:
8% 1. Decompose states into transient and recurrent classes via SCC
9% 2. For transient states: solve n * Q_tt = -p0_t for expected sojourn
10% 3. Compute hitting probabilities: h = n * Q_ta + p0_r
11% 4. For each recurrent class: solve pi_c * Q_cc = 0, scale by hitting prob
12%
13% Input:
14% Q: infinitesimal generator matrix
15% pi0: initial distribution vector, set to [] if not available
16% options: struct where options.tol sets the tolerance
17%
18% Output:
19% pi: limiting distribution (1 x N)
20% - For an ergodic CTMC, this is the unique limiting distribution.
21% - For a reducible CTMC:
22% - if there is a single transient SCC then this is the limiting
23% distribution when starting uniformly within it.
24% - otherwise pi is the weighted average of pis rows.
25% pis: limiting distribution given initialization in a single SCC (numSCC x N)
26% pi0: starting distribution for each row of pis (numSCC x N)
27% scc: mapping of each state of Q to its SCC index
28% isrec: element i is true if SCC i is recurrent
29
30N = size(Q, 1);
31
32if nargin < 3
33 options = struct('tol', 1e-12);
34end
35
36if nargin < 2
37 pin = [];
38else
39 pin = pi0;
40end
41
42% Ensure valid generator
43Q = ctmc_makeinfgen(Q);
44
45% Build adjacency from off-diagonal positive entries
46Adj = Q;
47Adj(1:N+1:end) = 0;
48[scc, isrec] = stronglyconncomp(Adj > 0);
49numSCC = max(scc);
50
51% Irreducible case: use standard solver
52if numSCC == 1
53 pi = ctmc_solve(Q, struct('force',true));
54 pis = pi;
55 pi0 = [];
56 return
57end
58
59% Build SCC index sets
60scc_idx = cell(1, numSCC);
61for i = 1:numSCC
62 scc_idx{i} = find(scc == i);
63end
64
65% Classify SCCs
66trans_scc_ids = find(~isrec);
67rec_scc_ids = find(isrec);
68
69% Gather ordered state indices
70trans_states = [];
71for i = trans_scc_ids
72 trans_states = [trans_states, scc_idx{i}];
73end
74trans_states = sort(trans_states);
75
76rec_states = [];
77for i = rec_scc_ids
78 rec_states = [rec_states, scc_idx{i}];
79end
80rec_states = sort(rec_states);
81
82nt = length(trans_states);
83nr = length(rec_states);
84
85% Extract Q sub-blocks (only when transient states exist)
86if nt > 0
87 Q_tt = Q(trans_states, trans_states);
88 Q_ta = Q(trans_states, rec_states);
89end
90
91% Compute per-SCC limiting distributions
92pis = zeros(numSCC, N);
93pi0 = zeros(numSCC, N);
94
95for s = 1:numSCC
96 % Starting distribution: uniform within SCC s
97 p0 = zeros(1, N);
98 p0(scc_idx{s}) = 1 / length(scc_idx{s});
99 pi0(s, :) = p0;
100
101 % Compute absorption probabilities into recurrent states
102 hit = zeros(1, nr);
103 if nt > 0
104 p0_t = p0(trans_states);
105 if any(abs(p0_t) > 0)
106 % Solve n * Q_tt = -p0_t for expected sojourn in transient states
107 % Q_tt is non-singular (Hurwitz) for transient states
108 hit = (-p0_t) / Q_tt * Q_ta;
109 end
110 end
111
112 % Add initial mass already in recurrent states
113 hit = hit + p0(rec_states);
114
115 % Solve steady-state per recurrent class, scaled by hitting probability
116 for c = rec_scc_ids
117 idx_c = scc_idx{c};
118 [~, loc] = ismember(idx_c, rec_states);
119 reachprob = sum(hit(loc));
120 if reachprob < 1e-15
121 continue
122 end
123 if length(idx_c) == 1
124 % Absorbing state: hitting probability IS the final probability
125 pis(s, idx_c) = reachprob;
126 else
127 % Solve pi_c * Q_cc = 0 within this recurrent class
128 pi_c = ctmc_solve(Q(idx_c, idx_c), struct('force',true));
129 pis(s, idx_c) = pi_c * reachprob;
130 end
131 end
132end
133
134% Compute initial SCC probabilities for weighted average
135if isempty(pin)
136 pinl = ones(1, numSCC);
137 % Zero out SCCs containing states with zero column sums (no incoming)
138 col_sums = sum(abs(Q), 1);
139 for j = find(col_sums < 1e-12)
140 pinl(scc(j)) = 0;
141 end
142 if sum(pinl) > 0
143 pinl = pinl / sum(pinl);
144 else
145 pinl = ones(1, numSCC) / numSCC;
146 end
147else
148 pinl = zeros(1, numSCC);
149 for i = 1:numSCC
150 pinl(i) = sum(pin(scc_idx{i}));
151 end
152end
153
154% Weighted average over starting SCCs
155pi = zeros(1, N);
156for i = 1:numSCC
157 if pinl(i) > 0
158 pi = pi + pis(i, :) * pinl(i);
159 end
160end
161
162% Special case: single transient SCC without explicit initial distribution
163if isscalar(trans_scc_ids) && isempty(pin)
164 pi = pis(trans_scc_ids, :);
165end
166
167% Normalize
168total = sum(pi);
169if total > 0
170 pi = pi / total;
171end
172
173end