LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
mamap2m_fit.m
1function MMAP = mamap2m_fit(M1, M2, M3, GAMMA, P, F, B, S, fbsWeights)
2% Fits a MAPH(2,m) or A MAMAP(2,m) that matches the provided
3% characteristics. Three characteristics of the marked process are matched
4% either exactly or approximately. Among these, the class probabilities are
5% always matched exactly. The remaining two characteristics, by default,
6% are the forward and backward moments, unless the underlying AMAP(2) is
7% degenerate. Different pairs of characteristics can be chosen by
8% specifying higher weights in the optional parameter fbsWeights.
9%
10% Input
11% - M1,M2,M3: moments of the inter-arrival times
12% - GAMMA: auto-correlation decay rate of the inter-arrival times
13% - P: class probabilities
14% - F: first-order forward moments
15% - B: first-order backward moments
16% - fbsWeights: the weight assigned to forward moments, backward moments
17% and class transition probabilities (default: [1, 1, 1])
18% Output
19% - MMAP: fitted MAMAP[m]
20
21if nargin == 2
22 % by default equal weights
23 % when weights are equal, F+B is preferred over F+S and B+S
24 fbsWeights = [1 1 1];
25end
26
27fbWeights = [fbsWeights(1), fbsWeights(2)];
28fsWeights = [fbsWeights(1), fbsWeights(3)];
29bsWeights = [fbsWeights(2), fbsWeights(3)];
30
31gammatol = 1e-4;
32degentol = 1e-8;
33
34if length(P) > 2
35 fprintf('Fitting MAMAP(2,m): fitting F+B because m > 2\n');
36 MMAP = mamap2m_fit_gamma_fb(M1, M2, M3, GAMMA, P, F, B);
37 return;
38end
39
40% TODO: if it is a Poisson process, convert to a second-order Poisson
41% in order to have more degrees of freedom for the marked process
42
43% TODO: if it is a PH-Renewal process and forward moments are preferred,
44% convert to non-canonical form
45
46if abs(GAMMA) < gammatol
47 fprintf('Fitting MAMAP(2,m): fitting MAPH because gamma = %e\n', GAMMA);
48 MMAP = maph2m_fit(M1, M2, M3, P, B);
49 return;
50end
51
52% fit underlying AMAP(2)
53[~,MAPS] = amap2_fit_gamma(M1, M2, M3, GAMMA);
54
55% fit marked Poisson process
56if length(MAPS) == 1 && size(MAPS{1}{1},1) == 1
57 fprintf('Fitting MAMAP(2,m): fitting marked Poisson because the underlying process has one state\n');
58 MAP = MAPS{1};
59 % number of classes
60 m = length(P);
61 % fit class probabilities
62 MMAP = cell(1,2+m);
63 MMAP{1} = MAP{1};
64 MMAP{2} = MAP{2};
65 for c = 1:m
66 MMAP{2+c} = MMAP{2} .* P(c);
67 end
68 return;
69end
70
71MMAPS = cell(1,length(MAPS));
72ERRORS = zeros(1,length(MAPS));
73for j = 1:length(MAPS)
74 h1 = -1/MAPS{j}{1}(1,1);
75 h2 = -1/MAPS{j}{1}(2,2);
76 r1 = h1*MAPS{j}{1}(1,2);
77 r2 = h2*MAPS{j}{2}(2,2);
78 % detect and handle degenerate cases
79 degen = 1;
80 if GAMMA > 0
81 if r1 < degentol || (1-r2) < degentol
82 error('Should not happen');
83 elseif abs(h2-h1*r2) < degentol
84 MMAPS{j} = mamap22_fit_fs_multiclass(MAPS{j}, P, F, S, [], fsWeights);
85 elseif abs(h1 - h2 + h2*r1) < degentol
86 MMAPS{j} = mamap22_fit_bs_multiclass(MAPS{j}, P, B, S, [], bsWeights);
87 elseif (1-r1) < degentol
88 % non-canonical APH(2)
89 % only forward can be fitted
90 MMAPS{j} = mamap22_fit_fs_multiclass(MAPS{j}, P, F, S, [], fsWeights);
91 elseif r2 < degentol
92 % canonical APH(2)
93 MMAPS{j} = maph2m_fit_multiclass(MAPS{j}, P, B);
94 else
95 degen = 0;
96 end
97 else
98 if (1-r2) < degentol
99 error('Should not happen');
100 elseif abs(h1 - h2 - h1*r1 + h1*r1*r2) < degentol
101 MMAPS{j} = mamap22_fit_fs_multiclass(MAPS{j}, P, F, S, [], fsWeights);
102 elseif abs(h1 - h2 + h2*r1) < degentol
103 MMAPS{j} = mamap22_fit_bs_multiclass(MAPS{j}, P, B, S, [], bsWeights);
104 elseif r2 < degentol && (1-r1) < degentol
105 % degenerate canonical APH(2)
106 MMAPS{j} = maph2m_fit_multiclass(MAPS{j}, P, B);
107 elseif r2 < degentol
108 if fbsWeights(1) >= fbsWeights(2)
109 % fit forward or sigma
110 MMAPS{j} = mamap22_fit_fs_multiclass(MAPS{j}, P, F, S, [], fsWeights);
111 else
112 % fit backward or sigma
113 MMAPS{j} = mamap22_fit_bs_multiclass(MAPS{j}, P, B, S, [], bsWeights);
114 end
115 else
116 degen = 0;
117 end
118 end
119 % handle non-degenerate cases according to user preference
120 if ~degen
121 if fbsWeights(1) >= fbsWeights(3) && fbsWeights(2) >= fbsWeights(3)
122 % prefer forward and backward
123 MMAPS{j} = mamap2m_fit_fb_multiclass(MAPS{j}, P, F, B, [], fbWeights);
124 elseif fbsWeights(1) >= fbsWeights(2)
125 % prefer forward and sigma
126 MMAPS{j} = mamap22_fit_fs_multiclass(MAPS{j}, P, F, S, [], fsWeights);
127 else
128 % prefer backward and sigma
129 MMAPS{j} = mamap22_fit_bs_multiclass(MAPS{j}, P, B, S, [], bsWeights);
130 end
131 end
132 % compute fitting error
133 fF = mmap_forward_moment(MMAPS{j}, 1);
134 fB = mmap_backward_moment(MMAPS{j}, 1);
135 fS = mmap_sigma(MMAPS{j});
136 ERRORS(j) = fbsWeights(1) * (F(1)/fF(1)-1)^2 + ...
137 fbsWeights(2) * (B(1)/fB(1)-1)^2 + ...
138 fbsWeights(3) * (S(1,1)/fS(1,1)-1)^2;
139end
140
141% pick the best fit
142[~,BEST] = min(ERRORS);
143MMAP = MMAPS{BEST};
144
145end