LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
QBD_IS.m
1function [G,R,U]=QBD_IS(A0,A1,A2,varargin)
2%QBD_IS Invariant Subspace for Quasi-Birth-Death Markov Chains [Akar, Sohraby]
3%
4% DISCRETE TIME CASE:
5%
6% G=QBD_IS(A0,A1,A2) computes the minimal nonnegative solution to the
7% matrix equation G = A0 + A1 G + A2 G^2, where A,B and C are square
8% nonnegative matrices, with (A0+A1+A2) irreducible and stochastic
9%
10% [G,R]=QBD_IS(A0,A1,A2) also provides the minimal nonnegative solution
11% to the matrix equation R = A2 + R A1 + R^2 A0
12%
13% [G,R,U]=QBD_IS(A0,A1,A2) also provides the minimal nonnegative solution
14% to the matrix equation U = A1 + A2 (I-U)^(-1) A0
15%
16% CONTINUOUS TIME CASE:
17%
18% G=QBD_IS(A0,A1,A2) computes the minimal nonnegative solution to the
19% matrix equation 0 = A0 + A1 G + A2 G^2, where A,B and C are square
20% nonnegative matrices, with (A0+A1+A2) having row sums equal to zero
21%
22% [G,R]=QBD_IS(A0,A1,A2) also provides the minimal nonnegative solution
23% to the matrix equation 0 = A2 + R A1 + R^2 A0
24%
25% [G,R,U]=QBD_IS(A0,A1,A2) also provides the minimal nonnegative solution
26% to the matrix equation U = A1 + A2 (-U)^(-1) A0
27%
28% Optional Parameters:
29%
30% MaxNumIt: Maximum number of iterations (default: 50)
31% Verbose: The residual error is printed at each step when set to 1,
32% (default:0)
33% Mode: 'MSignStandard' uses the matrix sign approach
34% 'MSignBalzer' uses the matrix sign approach with Balzer acceleration
35% 'Schur' relies on an ordered Schur decomposition to find the
36% invariant subspace (default: Schur)
37% RAPComp: set to 1 if the QBD has RAP components
38
39OptionNames=['Mode ';
40 'Verbose ';
41 'MaxNumIt ';
42 'RAPComp '];
43OptionTypes=['char ';
44 'numeric';
45 'numeric';
46 'numeric'];
47
48OptionValues{1}=['MSignStandard';
49 'MSignBalzer ';
50 'Schur '];
51
52options=[];
53for i=1:size(OptionNames,1)
54 options.(deblank(OptionNames(i,:)))=[];
55end
56
57% Default settings
58options.Mode='Schur';
59options.Verbose=0;
60options.MaxNumIt=50;
61options.RAPComp=0;
62
63% Parse Optional Parameters
64options=ParseOptPara(options,OptionNames,OptionTypes,OptionValues,varargin);
65
66if(~options.RAPComp)
67 % Convert to discrete time problem, if needed
68 m=size(A1,1);
69 continues=0;
70 if (sum(diag(A1)<0)) % continues time
71 continues=1;
72 lamb=max(-diag(A1));
73 A0=A0/lamb;
74 A1=A1/lamb+eye(m);
75 A2=A2/lamb;
76 end
77
78 % Parse Parameters
79 QBD_ParsePara(A0,A1,A2);
80else
81 % Parse Parameters
82 QBD_RAP_ParsePara(A0,A1,A2);
83
84 % Convert to discrete time problem - uniformization
85 m=size(A1,1);
86 continues=1;
87 lamb=max(-diag(A1));
88 A0=A0/lamb;
89 A1=A1/lamb+eye(m);
90 A2=A2/lamb;
91end
92
93% check whether G is known explicitly
94[G,R,U]=QBD_EG(A0,A1,A2,options.Verbose,nargout);
95if (~isempty(G))
96 return
97end
98
99epsilon=10^(-12);
100f=2;
101m=size(A1,1);
102
103
104% Step 1, F(z)=zD(z)-N(z), with A(z) = D^(-1)(z)N(z)
105% For QBD: A(z) is a polynomial matrix in z with degree f. Thus, D(z)=I.
106theta=stat(A0+A1+A2);
107drift=theta*sum(A0,2)-theta*sum(A2,2);
108F{1}=-A0;
109F{2}=eye(m)-A1;
110F{3}=-A2;
111
112% Step 2, H(s)=sum_{i=0}^f F_i (1-s)^(f-i) (1+s)^i = sum_{i=0}^f H_i s^i
113for i=0:f
114 H{i+1}=zeros(m,m);
115end
116for i=0:f
117 temp1=[1 -1];
118 temp2=[1 1];
119 con1=1;
120 con2=1;
121 for j=1:f-i
122 con1=conv(con1,temp1);
123 end
124 for j=1:i
125 con2=conv(con2,temp2);
126 end
127 contrib=conv(con1,con2);
128 for j=0:f
129 H{j+1}=H{j+1}+contrib(j+1)*F{i+1};
130 end
131end
132
133clear F;
134
135% Step 3, \hat{H}_i = H_f^-1*H_i
136H{f+1}=inv(H{f+1});
137for i=0:f-1
138 hatH{i+1}=H{f+1}*H{i+1};
139end
140
141% Step 4, y, xT
142y=[ones(m,1); zeros(m*(f-1),1)];
143x0T=[zeros(1,m) 1] / [hatH{1} ones(m,1)];
144for i=1:f-1
145 xT(1,(i-1)*m+1:i*m)=x0T*hatH{i+1};
146end
147xT(1,(f-1)*m+1:f*m)=x0T;
148
149% Step 5, E_m in Zold
150Zold=zeros(m*f,m*f);
151for i=1:f-1
152 Zold((i-1)*m+1:i*m,i*m+1:(i+1)*m)=eye(m);
153end
154for i=0:f-1
155 Zold(m*(f-1)+1:m*f,i*m+1:(i+1)*m)=-hatH{i+1};
156end
157y=y/(xT*y);
158Zold=Zold-sign(drift)*y*xT;
159
160if ( exist('ordschur') ~= 5 || ~strcmp(options.Mode,'Schur'))
161 % Step 6, classic matrix sign function algorithm
162 if (strcmp('Schur',options.Mode))
163 fprintf('Ordschur not supported by current MATLAB version\n');
164 fprintf('An automatic switch is performed to the MSignBalzer Mode\n');
165 drawnow;
166 end
167 numit=0;
168 check=1;
169 while (check > epsilon && numit < options.MaxNumIt)
170 numit=numit+1;
171 if (strcmp(options.Mode,'MSignStandard'))
172 determ=1/2;
173 else % options.Mode = 'MSignBalzer' or switched from 'Schur' Mode
174 determ=(1+abs(det(Zold))^(1/(m*f)))^(-1);
175 determ=min(determ,1-10^(-3));
176 end
177 Znew=determ*Zold+(1-determ)*inv(Zold);
178 check=norm(Znew-Zold,1)/norm(Zold,1);
179 if (options.Verbose==1)
180 fprintf('Check after %d iterations: %d\n',numit,check);
181 drawnow;
182 end
183 Zold=Znew;
184 end
185 if (numit == options.MaxNumIt && check > epsilon)
186 warning('Maximum Number of Iterations %d reached: T may not have m columns',numit);
187 end
188 % Step 7,
189 T=orth(Znew-eye(m*f));
190else
191 % Theorem 7 (Akar, Oguz, Sohraby -- SM 2000)
192 [T,D]=schur(Zold);
193 [T,~]=ordschur(T,D,'lhp');
194 T=T(:,1:m);
195end
196
197% Step 8,
198G=(T(1:m,:)+T(m+1:2*m,:))*inv(T(1:m,:)-T(m+1:2*m,:));
199
200if (options.Verbose==1)
201 res_norm=norm(G-A0-(A1+A2*G)*G,inf);
202 fprintf('Final Residual Error for G: %d\n',res_norm);
203end
204
205% Compute R
206if (nargout > 1)
207 R=A2*(eye(m)-(A1+A2*G))^(-1);
208 if (options.Verbose==1)
209 res_norm=norm(R-A2-R*(A1+R*A0),inf);
210 fprintf('Final Residual Error for R: %d\n',res_norm);
211 end
212end
213
214% Compute U
215if (nargout > 2)
216 U=A1+R*A0;
217 if (options.Verbose==1)
218 res_norm=norm(U-A1-A2*(eye(m)-U)^(-1)*A0,inf);
219 fprintf('Final Residual Error for U: %d\n',res_norm);
220 end
221 if (continues)
222 U=lamb*(U-eye(m));
223 end
224end
225