LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
kpcfit_ph_auto.m
1function [PH] = kpcfit_ph_auto(E,options)
2% [PH] = kpcfit_ph_auto(E, options)
3%
4% DESCRIPTION
5% Automatic fit a phase-type (PH) distribution using exact or approximate
6% moment matching. Approximations are based on the Kronecker Product
7% Composition method (KPC).
8%
9% INPUT
10% E - vector of moments, E(k) = E[X^k]
11% options - ph fitting
12%
13% OUTPUT
14% PH{i,1}: the i-th fitted PH described as a PH-renewal process specified
15% in the (D0,D1) notation. The (alpha,T) representa
16%
17% EXAMPLE
18%
19% S; % trace
20% for k=1:16 E(k) = mean(S.^k); end
21% options = kpcfit_ph_options(E,'MinNumStates',2,'MaxNumStates',8);
22% [PH]=kpcfit_ph_auto(E,options)
23% kpcfit_ph_summary
24% i=1; [alpha,T]=map2ph(PH{i,1}) % (alpha,T) distribution of i-th fitted PH
25
26if options.Verbose
27 fprintf('kpcfit_ph: version %s\n',kpcfit_version())
28 fprintf('kpcfit_ph: type "help kpcfit_ph_options" for options\n')
29end
30
31%% initialization
32% scale moments such that E[X]=1
33Eunscaled = E;
34Escale = [];
35for k = 1:length(E)
36 Escale(k) = E(1).^k;
37end
38E = (E(:) ./ Escale(:))';
39
40%% determine best exact fit -- Prony's method
41if options.Verbose
42 fprintf('\nkpcfit_ph: starting exact fitting methods\n')
43end
44PH_EXACT = kpcfit_ph_exact(E, options);
45
46%% determine best approximate fit -- moment space
47if options.Verbose
48 fprintf('\nkpcfit_ph: starting approximate fitting method -- moment space (KPC-based)\n')
49 fprintf('\t\t\tRUN\t\tORD\tDIST\t\tRESULT')
50end
51PH_APX_MS = {};
52
53tic;
54for J = max([2,ceil(log2(options.MinNumStates))]):ceil(log2(options.MaxNumStates)) % for all number of MAPs to compose
55 % initialize
56 mindist = Inf;
57 best = 0;
58 bestidx = 0;
59
60 % random runs
61 randomruns = ceil(2*options.Runs/3); % number of runs with random initial point
62 for run = 1:options.Runs
63 if options.Verbose
64 if run > randomruns
65 fprintf(1,'\n\t\t\t%dr\t',run); %run number
66 else
67 fprintf(1,'\n\t\t\t%d\t',run); %run number
68 end
69 fprintf(1,'\t%d',2^J); % number of states
70 end
71 %try
72 if run > randomruns && best > 0
73 x0 = PH_APX_MS{bestidx,3}; %x0 = x0 * mean(x0)*1e-2*rand(size(x0));
74 else
75 x0 = []; % automatic
76 end
77 [APX, ~, x] = kpcfit_ph_search(E,J,options,x0);
78
79 if map_isfeasible(APX)
80 PH_APX_MS = store_result(PH_APX_MS,APX);
81 else
82 if options.Verbose
83 fprintf(1,'\t\t\t\tx infeasible result.');
84 end
85 end
86 %catch
87 % fprintf(1,'\t\t\t\tx optimization solver failed.');
88 %end
89 end
90end
91
92%% determine best approximate fit -- parameter space
93
94SCV = (E(2)-E(1)^2)/E(1)^2;
95if SCV >= 1
96 if options.Verbose
97 fprintf('\n\nkpcfit_ph: starting approximate fitting method -- hyper-exponential parameter space (non-KPC-based)\n')
98 end
99 PH_APX_PS = {};
100 if options.Verbose
101 fprintf('\t\t\tRUN\t\tORD\tDIST\t\tRESULT')
102 end
103 for J = max([2,ceil(log2(options.MinNumStates))]):ceil(log2(options.MaxNumStates)) % for all number of MAPs to compose
104 % initialize
105 mindist = Inf;
106 best = 0;
107 bestidx = 0;
108
109 % random runs
110 randomruns = ceil(2*options.Runs/3); % number of runs with random initial point
111 for run = 1:options.Runs
112 if options.Verbose
113 fprintf(1,'\n\t\t\t%d\t',run); %run number
114 fprintf(1,'\t%d',2^J); % number of states
115 end
116 try
117 [APX,~,x] = psfit_hyperexp_weighted(E,2^J);
118
119 if map_isfeasible(APX)
120 PH_APX_PS = store_result(PH_APX_PS,APX);
121 else
122 if options.Verbose
123 fprintf(1,'\t\t\t\tinfeasible result.');
124 end
125 end
126 catch
127 if options.Verbose
128 fprintf(1,'\t\t\t\tx optimization solver failed.');
129 end
130 end
131 end
132 end
133else
134 fprintf('\n\nkpcfit_ph: starting approximate fitting method -- aph parameter space (non-KPC-based)\n')
135 PH_APX_PS = {};
136 fprintf('\t\t\tRUN\t\tORD\tDIST\t\tRESULT')
137 for J = max([2,ceil(log2(options.MinNumStates))]):ceil(log2(options.MaxNumStates)) % for all number of MAPs to compose
138 % initialize
139 mindist = Inf;
140 best = 0;
141 bestidx = 0;
142
143 % random runs
144 randomruns = ceil(2*options.Runs/3); % number of runs with random initial point
145 for run = 1:options.Runs
146 fprintf(1,'\n\t\t\t%d\t',run); %run number
147 fprintf(1,'\t%d',2^J); % number of states
148 try
149 [APX,~,x] = psfit_aph_weighted(E,2^J);
150
151 if map_isfeasible(APX)
152 PH_APX_PS = store_result(PH_APX_PS,APX);
153 else
154 fprintf(1,'\t\t\t\tinfeasible result.');
155 end
156 catch
157 if options.Verbose
158 fprintf(1,'\t\t\t\tx optimization solver failed.');
159 end
160 end
161 end
162 end
163end
164
165%% ensure that all returned results match at least the mean
166nExactResults = size(PH_EXACT,1);
167nApproxMS = size(PH_APX_MS,1);
168%nApproxPS = size(PH_APX_PS,1);
169
170for j = 1:size(PH_EXACT,1)
171 PH{j,1} = map_scale(PH_EXACT{j},Eunscaled(1));
172 PH{j,2} = dist_fun(E,map_moment(PH_EXACT{j},1:length(E)));
173 PH{j,3} = [];
174 PH{j,4} = 'exact';
175end
176for j = 1:size(PH_APX_MS,1)
177 PH{nExactResults + j,1} = map_scale(PH_APX_MS{j,1},Eunscaled(1));
178 PH{nExactResults + j,2} = dist_fun(E,map_moment(PH_APX_MS{j,1},1:length(E)));
179 PH{nExactResults + j,3} = PH_APX_MS{j,3};
180 PH{nExactResults + j,4} = 'approx_moment_space';
181end
182for j = 1:size(PH_APX_PS,1)
183 PH{nExactResults + nApproxMS + j,1} = map_scale(PH_APX_PS{j,1},Eunscaled(1));
184 PH{nExactResults + nApproxMS + j,2} = dist_fun(E,map_moment(PH_APX_PS{j,1},1:length(E)));
185 PH{nExactResults + nApproxMS + j,3} = PH_APX_PS{j,3};
186 PH{nExactResults + nApproxMS + j,4} = 'approx_param_space';
187end
188
189if options.Verbose
190 fprintf(sprintf('\n\nReturned %d PH distribution(s) by exact moment matching.',nExactResults));
191 fprintf(sprintf('\nReturned %d PH distribution(s) by approximate moment matching (moment space).',size(PH_APX_MS,1)));
192 fprintf(sprintf('\nReturned %d PH distribution(s) by approximate moment matching (parameter space).',size(PH_APX_PS,1)));
193 fprintf('\n\n');
194end
195return
196
197 function dist = dist_fun(Eref,Eapx)
198 % weight function
199 w = (log10(Eref(end)).^(1/length(Eref))).^-(1:length(Eref));
200 dist = w*abs( log10(Eref) - log10(Eapx) )'; %%does not work!
201 end
202
203 function PH_APX_MS = store_result(PH_APX_MS,APX)
204 dist = dist_fun(E, map_moment(APX,1:length(E)));
205 if isnan(dist)
206 if options.Verbose
207 fprintf(1,'\t-------- \tx no solution found.',dist);
208 end
209 else
210 if options.Verbose
211 fprintf(1,'\t%6.6f',dist);
212 end
213 if min([Inf,dist]) < mindist
214 % fprintf(1,'+',length(PH_APX_MS));
215 if best == 0
216 PH_APX_MS{end+1,1} = APX;
217 if length(APX{1}) > 2^J
218 if options.Verbose
219 fprintf(1,sprintf('\to PH(%d) saved. \t elapsed: %0.3fs [order automatically increased to match second moment]',length(APX{1}),toc));
220 end
221 else
222 if options.Verbose
223 fprintf(1,sprintf('\to PH(%d) saved. \t elapsed: %0.3fs',length(APX{1}),toc));
224 end
225 end
226 else
227 PH_APX_MS{end,1} = APX;
228 if options.Verbose
229 fprintf(1,sprintf('\t+ PH(%d) updated. \t elapsed: %0.3fs',length(APX{1}),toc));
230 end
231 end
232 best = run;
233 bestidx = size(PH_APX_MS,1);
234 mindist = min([Inf,dist]);
235 PH_APX_MS{end,2} = dist;
236 PH_APX_MS{end,3} = x;
237 else
238 if options.Verbose
239 fprintf(1,sprintf('\t\t\t\t\t\t elapsed: %0.3fs',toc));
240 end
241 end
242 end
243 end
244end
245
246