1function [mmap,fF,fS,exact] = mamap22_fit_fs_multiclass(
map,p,F,S,classWeights,fsWeights,adjust)
2% Performs approximate fitting of a
MMAP given the underlying AMAP(2),
3% the
class probabilities (always fitted exactly), the forward moments,
4% and the one-step
class transition probabilities.
6% -
map: second-order AMAP underlying the MAMAP[2]
7% - p: vector of
class probabilities
8% - F: vector of forward moments
9% - S: matrix of the
class transition probabilities
10% - classWeights: optional vector of weights
for each class
11% - fsWeights: optional 2-vector of weights of forward moments and the
12%
class transition probabilities
14% - mmap: fitted MAMAP[2]
15% - fF: vector of optimal feasible forward moments
16% - fS: vector of optimal feasible
class transition probabilities
18if (size(
map{1},1) ~= 2)
19 error('Underlying MAP must be of second-order.');
22 error('Underlying MAP must be acyclic');
26elseif (
map{2}(1,1) == 0)
29 error('Underlying MAP must be in canonical acyclic form');
32fprintf('Fitting MAMAP(2,2) F+S: form = %d\n', form);
37 error('Fitting MAMAP(2,2) F+S: only two
classes supported');
40% default weights to use in the objective function
41if nargin < 5 || isempty(classWeights)
42 classWeights = ones(k,1);
44if nargin < 6 || isempty(fsWeights)
45 fsWeights = ones(2,1);
56% get parameters of the underlying AMAP(2)
62% set tolerance constants
67% generate required coefficients
69 [G,U,Y] = mamap2m_can1_coefficients(h1,h2,r1,r2);
72 [E,V,Z] = mamap2m_can2_coefficients(h1,h2,r1,r2);
78if (form == 1 && (r1 < degentol || r2 > 1-degentol || abs(h1 - h2 + h2*r1) < degentol )) || ...
79 (form == 2 && (r2 > 1-degentol || abs(h1 - h2 + h2*r1) < degentol ))
81 % Transform into a valid second-order Poisson process by perturbing
82 % degenerate parameters to maintain 2-state structure
83 fprintf(
'Fitting MAMAP(2,2) F+S: perturbing degenerate Poisson to second-order\n');
90 if abs(h1 - h2 + h2*r1) < degentol
91 h1 = h2 * (1 - r1) + degentol;
93 % Reconstruct MAP with perturbed parameters
95 map{1} = [-1/h1, r1/h1; 0, -1/h2];
96 map{2} = [(1-r1)/h1, 0; r2/h2, (1-r2)/h2];
98 mmap{1} =
map{1}; mmap{2} =
map{2};
99 [G,U,Y] = mamap2m_can1_coefficients(h1,h2,r1,r2);
102 map{1} = [-1/h1, r1/h1; 0, -1/h2];
103 map{2} = [0, (1-r1)/h1; r2/h2, (1-r2)/h2];
105 mmap{1} =
map{1}; mmap{2} =
map{2};
106 [E,V,Z] = mamap2m_can2_coefficients(h1,h2,r1,r2);
110elseif form == 2 && r2 < degentol && abs(1-r1) < degentol
112 % DEGENERATE PHASE_TYPE
113 fprintf(
'Fitting MAMAP(2,2) F+S: detected degenerate phase-type form\n');
115 % only one degree of freedom: match
class probabilites
120elseif form == 1 && r2 < degentol
122 % CANONICAL PHASE_TYPE
123 fprintf(
'Fitting MAMAP(2,2) F+S: detected canonical phase-type form\n');
125 % convert to phase-type
130 % detect hypoexponential (SCV <= 1) and convert to non-canonical form
132 if scv < 1.0 + degentol
133 fprintf('Fitting MAMAP(2,2) F+S: hypoexponential detected (SCV=%.4f), converting to non-canonical form\n', scv);
137 % the forward moments are always equal to the ordinary moments
138 % since the backward moments are not provided as arguments,
139 % we set the backward moments to the ordinary moments as well
141 B(1) = map_mean(
aph);
142 B(2) = map_mean(
aph);
143 warning('Fitting MAMAP(2,2) F+S: setting backward moments to ordinary moments, you should try to fit B+S');
145 mmap = maph2m_fit_multiclass(
aph, p, B, classWeights);
147 fF = mmap_forward_moment(mmap, 1);
148 fS = mmap_sigma(mmap);
152elseif (form == 1 && abs(1-r1) < degentol) || ...
153 (form == 2 && abs(1-r1) < degentol)
155 % NON-CANONICAL PHASE_TYPE
156 fprintf('Fitting MAMAP(2,2) F+S: detected non-canonical phase-type form\n');
158 % two degrees of freedom: fit p and F, not S
159 [q1,q2,q3] = fit_can1_degen(F(1));
161 % if unfeasible, perform approximate fitting
162 if ~(isfeasible(q1) && isfeasible(q2) && isfeasible(q3))
163 % four inequalities, one variable
167 q2_F = - p(1) / (h1 * (r2-1));
168 q2_0 = + p(1) * h2 / (h1 * (r2-1));
169 q3_F = - p(1) / (h1 * r2);
170 q3_0 = + p(1) * (h1 + h2) / (h1 * r2);
187 fF1 = solve_quadprog();
188 % compute coefficient
189 [q1, q2, q3] = fit_can1_degen(fF1);
192elseif form == 2 && r2 < degentol
194 % DEGENERATE CASE FOR gamma < 0
195 fprintf('Fitting MAMAP(2,2) F+S: detected degenerate
MMAP form\n');
197 if fsWeights(1) > fsWeights(2)
198 fprintf('Fitting MAMAP(2,2) F+S: fitting F\n');
199 [q1,q2,q3] = fit_can2_degen_forward(F(1));
200 if ~(isfeasible(q1) && isfeasible(q2) && isfeasible(q3))
201 % four inequalities, one variable
205 q1_F = - p(1) * (r1 - 2) / ((r1-1)*(h1-h2+h2*r1));
206 q1_0 = + p(1) * (r1 - 2)*(h1 + h2*r1) / ((r1-1)*(h1-h2+h2*r1));
207 q2_F = - p(1) * (r1 - 2) / ((h1-h2+h2*r1));
208 q2_0 = + p(1) * (r1 - 2)*h2 / ((h1-h2+h2*r1));
225 fF1 = solve_quadprog();
226 % compute coefficient
227 [q1,q2,q3] = fit_can2_degen_forward(fF1);
230 fprintf('Fitting MAMAP(2,2) F+S: fitting S\n');
231 [q1,q2,q3] = fit_can2_degen_transition(S(1,1));
232 if ~(isfeasible(q1) && isfeasible(q2) && isfeasible(q3))
233 % nonlinear constraints: use yalmip
234 varS11 = sdpvar(1,1);
235 options = sdpsettings('verbose',0);
236 boundS11 = [varS11 >= 0, varS11 <= p(1)^2];
237 q1lb = varS11 >= p(1)^2 * (1 - (1-r1)^2);
238 q2ub = varS11 >= p(1)^2 - (1 - p(1))^2;
239 constraints = [boundS11, q1lb, q2ub];
240 obj = (varS11 - S(1,1))^2;
241 sol = solvesdp(constraints, obj, options);
243 error('Fitting MAMAP(2,2) F+S: YALMIP = %s\n', sol.info);
245 fprintf('Fitting MAMAP(2,2) F+S: YALMIP objective = %f\n',
double(obj));
247 % compute coefficient
248 [q1,q2,q3] = fit_can2_degen_transition(
double(varS11));
254 % FULL FORM or "GOOD" poisson process
256 if (form == 1 && abs(h2-h1*r2) < degentol) || (form == 2 && abs(h1 - h2 - h1*r1 + h1*r1*r2) < degentol)
257 fprintf('Fitting MAMAP(2,2) F+S: detected fittable Poisson process\n')
260 [q1,q2,q3] = fit(F(1), S(1,1));
262 % options used to solve the nonlinear problem
263 use_linprog_prepass = 1;
265 optim_space = 'hybrid';
266 %optim_space = 'hybrid-z';
267 %optim_space = 'param';
268 %optim_space = '
char';
270 yalmip_nonlinear_opt = sdpsettings(...
272 'solver','bmibnb',...
276 'bmibnb.relgaptol',1e-3,...
277 'bmibnb.absgaptol',1e-6,...
278 'bmibnb.maxiter',1000,...
279 'bmibnb.lowrank', 1,...
280 'bmibnb.lpreduce',1,... % without this, crappy lower bounds for some problems
281 'bmibnb.pdtol',-1e-8,... % x >= if x > -pdtol
282 'bmibnb.eqtol',+1e-10,... % x == 0 if abs(x) < +eqtol
283 'bmibnb.roottight',1,...
284 'bmibnb.uppersolver','fmincon-standard',...
285 'fmincon.Algorithm','sqp',...
286 'fmincon.TolCon',1e-6,...
287 'fmincon.TolX',1e-6,...
288 'fmincon.GradObj','on',...
289 'fmincon.GradConstr','on',...
290 'fmincon.Hessian','off',...
291 'fmincon.LargeScale','off');
296 if isfeasible(q1) && isfeasible(q2) && isfeasible(q3)
297 fprintf('Fitting MAMAP(2,2) F+S: exact fit found\n');
300 [q1,q2,q3] = solve_nonlinear(optim_space);
303end % end if/then/else for degenerate forms
305% check parameter feasibility
306if adjust && ~(isfeasible(q1) && isfeasible(q2) && isfeasible(q3))
307 error('Fitting MAMAP(2,2) F+S: Feasibility could not be restored: q1 = %e, q2 = %e, q3 = %e', q1, q2, q3);
309% parameters feasible within feastol: restrict to [0,1]
314% compute D11,D12,...,D1k
316 mmap{2+1} = mmap{2} .* [q1 0; q2 q3];
317 mmap{2+2} = mmap{2} .* [1-q1 0; 1-q2 1-q3];
319 mmap{2+1} = mmap{2} .* [0 q1; q2 q3];
320 mmap{2+2} = mmap{2} .* [0 1-q1; 1-q2 1-q3];
323fF = mmap_forward_moment(mmap, 1);
324fS = mmap_sigma(mmap);
326 function feas = isfeasible(q)
327 feas = q >= -feastol && q <= (1+feastol);
330 function qfix = fix(q)
331 qfix = max(min(q,1),0);
334 function [q1,q2,q3] = fit_can1(vF1,vS11)
335 denum = p(1) * (U(5) * vF1 + U(6));
336 if abs(denum) < denumtol
341 q2 = (U(1) * vF1^2 * p(1)^2 + U(2) * vF1 * p(1)^2 + U(3) * vS11 + U(4) * p(1)^2)/denum;
342 q1 = -(G(15)*p(1) - vF1*G(3)*p(1) + (G(3)*G(14) - G(2)*G(15)) * q2)/Y(2);
343 q3 = +(G(13)*p(1) - vF1*G(1)*p(1) + (G(1)*G(14) - G(2)*G(13)) * q2)/Y(2);
347 function [q1,q2,q3] = fit_can1_degen(vF1)
348 q1 = p(1); % meangingless
if really degenerate
349 q2 = p(1) * (h2 - vF1) / (h1 * (r2-1));
350 q3 = p(1) * (h1 + h2 - vF1) / (h1 * r2);
353 function [q1,q2,q3] = fit_can2_degen_forward(vF1)
354 q1 = + p(1) * (r1 - 2)*(h1 + h2*r1 - vF1)/((r1-1)*(h1-h2+h2*r1));
355 q2 = - p(1) * (vF1 - h2)*(r1 - 2)/((h1-h2+h2*r1));
356 q3 = p(1); % meangingless
if really degenerate
359 function [q1,q2,q3] = fit_can2_degen_transition(vS11)
360 q1 = p(1) + (p(1)^2 - vS11)^(1/2)/(r1 - 1);
361 q2 = p(1) + (p(1)^2 - vS11)^(1/2);
362 q3 = p(1); % meangingless
if really degenerate
365 function [q1,q2,q3] = fit_can2(vF1,vS11)
366 denum = (V(5) * vF1 * p(1) + V(6) * p(1));
367 if abs(denum) < denumtol
372 q3 = (V(1) * vF1^2 * p(1)^2 + V(2) * vF1 * p(1)^2 + V(3) * p(1)^2 + V(4) * vS11)/denum;
373 q1 = -(E(13)*p(1) - vF1*E(2)*p(1) + (E(2)*E(14) - E(3)*E(13)) * q3)/Z(2);
374 q2 = +(E(12)*p(1) - vF1*E(1)*p(1) + (E(1)*E(14) - E(3)*E(12)) * q3)/Z(2);
378 function [pexp,pFexp,Sexp] = get_characteristics(q1,q2,q3)
380 pexp = G(1)*q1+G(2)*q2+G(3)*q3;
381 pFexp = G(13)*q1+G(14)*q2+G(15)*q3;
382 Sexp = G(4)*q1^2 + G(5)*q1*q2 + G(6)*q1*q3 + G(7)*q2^2 + G(8)*q2*q3 + G(9)*q3^2;
384 pexp = E(1)*q1+E(2)*q2+E(3)*q3;
385 pFexp = (E(12)*q1+E(13)*q2+E(14)*q3);
386 Sexp = E(4)*q1*q2 + E(5)*q1*q3 + E(6)*q2^2 + E(7)*q2*q3 + E(8)*q3^2;
390 function opt = get_yalmip_linear_options()
391 % we
do not set solver so that YALMIP picks the best
392 opt = sdpsettings(
'verbose',0);
395 function [q1,q2,q3,obj] = solve_nonlinear(optim_space)
397 fprintf(
'Fitting MAMAP(2,2) F+S: running approximate fitting (nonlinear)\n');
399 if strcmp(optim_space,
'param')
400 [q1,q2,q3,obj] = solve_nonlinear_param();
402 if strcmp(optim_space,
'hybrid')
403 solve_nonlinear_func = @solve_nonlinear_hybrid;
404 elseif strcmp(optim_space, 'hybrid-z')
405 solve_nonlinear_func = @solve_nonlinear_hybrid_z;
406 elseif strcmp(optim_space, '
char')
407 solve_nonlinear_func = @solve_nonlinear_char;
409 error('Fitting MAMAP(2,2) F+S: invalid value for option "optim_space"');
412 eobj = make_objective(M1, p(1)^2);
413 fprintf('Fitting MAMAP(2,2) F+S: F1 == M1, objective = %e\n', eobj);
415 [lfeas,lobj,lF1,lS11] = solve_nonlinear_func('left');
417 fprintf('Fitting MAMAP(2,2) F+S: F1 < M1, objective = %e\n', lobj);
419 fprintf('Fitting MAMAP(2,2) F+S: F1 < M1 infeasible\n');
422 [rfeas,robj,rF1,rS11] = solve_nonlinear_func('right');
424 fprintf('Fitting MAMAP(2,2) F+S: F1 > M1, objective = %e\n', robj);
426 fprintf('Fitting MAMAP(2,2) F+S: F1 > M1 infeasible\n');
429 if (~lfeas && ~rfeas) || (eobj < lobj && eobj < robj)
434 elseif ~rfeas || lobj < robj
436 [q1,q2,q3] = fit(lF1, lS11);
439 [q1,q2,q3] = fit(rF1, rS11);
445% used for non-degenerate cases: parmeter space
446 function [q1,q2,q3,obj] = solve_nonlinear_param()
447 % find initial feasible solution
451 [pexp,pFexp,Sexp] = get_characteristics(vq(1),vq(2),vq(3));
453 cstr_exact_p = pexp == p(1);
454 cstr_tighten_F = -z(1) <= (pFexp - p(1)*F(1))/(p(1)*F(1)) <= z(1);
455 cstr_tighten_S = -z(2) <= (Sexp - S(1,1))/S(1,1) <= z(2);
456 cstr = [0 <= vq(1) <= 1, ...
462 z(1) >= 0, z(2) >= 0];
464 zobj = fsWeights(1) * z(1)^2 + fsWeights(2) * z(2)^2;
466 zsol = solvesdp(cstr, zobj, yalmip_nonlinear_opt);
473 fprintf('Fitting MAMAP(2,2) F+S: solver (parameter space) objective = %e\n', obj);
476 save(fname,'
map','p','F','S');
477 error('Fitting MAMAP(2,2) F+S: solver (parameter space) error: %s, input saved to %s\n', zsol.info, fname);
481% used for non-degenerate cases: hybrid space
482 function [feas,bobj,bF1,bS11] = solve_nonlinear_hybrid(side)
488 q1exp = -(G(15)*p(1) - vF1*G(3)*p(1) + vq2*(G(14)*G(3)-G(15)*G(2)))/Y(2);
489 q2num = U(1)*vF1^2*p(1)^2 + U(2)*vF1*p(1)^2 + U(3)*vS11 + U(4)*p(1)^2;
490 q2den = p(1)*(U(5)*vF1 + U(6));
491 q3exp = +(G(13)*p(1) - vF1*G(1)*p(1) + vq2*(G(14)*G(1)-G(13)*G(2)))/Y(2);
492 hcstr = [0 <= q1exp <= 1, 0 <= vq2 <= 1, vq2*q2den/scale == q2num/scale, 0 <= q3exp <= 1, 0 <= vS11 <= 1, vF1 >= 0];
496 q1exp = -(E(13)*p(1) - vF1*E(2)*p(1) + vq3*(E(14)*E(2)-E(13)*E(3)))/Z(2);
497 q2exp = +(E(12)*p(1) - vF1*E(1)*p(1) + vq3*(E(14)*E(1)-E(12)*E(3)))/Z(2);
498 q3num = V(1)*vF1^2*p(1)^2 + V(2)*vF1*p(1)^2 + V(3)*p(1)^2 + V(4)*vS11;
499 q3den = p(1)*(V(5)*vF1 + V(6));
500 hcstr = [0 <= q1exp <= 1, 0 <= q2exp <= 1, vq3*q3den/scale == q3num/scale, 0 <= vq3 <= 1, 0 <= vS11 <= 1, vF1 >= 0];
502 if strcmp(side,'left')
503 hcstr = [hcstr, vF1 <= M1-tol];
505 hcstr = [hcstr, vF1 >= M1+tol];
507 hobj = make_objective(vF1,vS11);
508 hsol = solvesdp(hcstr, hobj, yalmip_nonlinear_opt);
510 fprintf('Fitting MAMAP(2,2) F+S: solver (hybrid space) objective = %f\n',
double(hobj));
515 elseif hsol.problem == 1 || hsol.problem == 12
516 fprintf('Fitting MAMAP(2,2) F+S: program
is infeasible\n');
523 save(fname,'
map','p','F','S');
524 error('Fitting MAMAP(2,2) F+S: solver (hybrid space) error: %s, input saved to %s\n', hsol.info, fname);
528% used for non-degenerate cases: hybrid space with z variable and implicit vS11
529 function [feas,bobj,bF1,bS11] = solve_nonlinear_hybrid_z(side)
534 q1exp = -(G(15)*p(1) - vF1*G(3)*p(1) + vq2*(G(14)*G(3)-G(15)*G(2)))/Y(2);
535 q2den = p(1)*(U(5)*vF1 + U(6));
536 q3exp = +(G(13)*p(1) - vF1*G(1)*p(1) + vq2*(G(14)*G(1)-G(13)*G(2)))/Y(2);
537 S11exp = (vq2*q2den - p(1)^2*(U(1)*vF1^2 + U(2)*vF1 + U(4)))/U(3);
538 tighten = -z <= S11exp - S(1,1) <= +z;
539 hcstr = [0 <= q1exp <= 1, 0 <= vq2 <= 1, 0 <= q3exp <= 1, vF1 >= 0, tighten];
542 q1exp = -(E(13)*p(1) - vF1*E(2)*p(1) + vq3*(E(14)*E(2)-E(13)*E(3)))/Z(2);
543 q2exp = +(E(12)*p(1) - vF1*E(1)*p(1) + vq3*(E(14)*E(1)-E(12)*E(3)))/Z(2);
544 q3den = p(1)*(V(5)*vF1 + V(6));
545 S11exp = (vq3*q3den - p(1)^2*(V(1)*vF1^2 + V(2)*vF1 + V(3)))/V(4);
546 tighten = -z <= S11exp - S(1,1) <= +z;
547 hcstr = [0 <= q1exp <= 1, 0 <= q2exp <= 1, 0 <= vq3 <= 1, vF1 >= 0, tighten];
549 if strcmp(side,'left')
550 hcstr = [hcstr, vF1 <= M1-tol];
552 hcstr = [hcstr, vF1 >= M1+tol];
554 hobj = fsWeights(1)*(vF1/F(1)-1)^2 + fsWeights(2)*(z/S(1,1))^2;
555 hsol = solvesdp(hcstr, hobj, yalmip_nonlinear_opt);
557 fprintf('Fitting MAMAP(2,2) F+S: solver (hybrid-z space) objective = %f\n',
double(hobj));
561 bS11 =
double(S11exp);
562 elseif hsol.problem == 1 || hsol.problem == 12
563 fprintf('Fitting MAMAP(2,2) F+S: program
is infeasible\n');
570 save(fname,'
map','p','F','S');
571 error('Fitting MAMAP(2,2) F+S: solver (hybrid-z space) error: %s, input saved to %s\n', hsol.info, fname);
575% used for non-degenerate cases: characteristic space
576 function [feas,bobj,bF1,bS11] = solve_nonlinear_char(side)
578 if strcmp(side,'left')
579 fprintf('Fitting MAMAP(2,2) F+S: optimizing for F11 < M1\n');
582 fprintf('Fitting MAMAP(2,2) F+S: optimizing for F11 > M1\n');
584 % find initial feasible solution
585 if use_linprog_prepass
586 % parameter variables vq and auxiliary variables
590 [pexp,pFexp,Sexp] = get_characteristics(vq(1),vq(2),vq(3));
593 cstr_side = pFexp <= pexp*(M1-tol);
595 cstr_side = pFexp >= pexp*(M1+tol);
597 cstr_exact_p = pexp == p(1);
598 cstr_tighten = -z <= pFexp-p(1)*F(1) <=z;
599 fcstr = [0 <= vq(1) <= 1, ...
606 % run linear programming solver
607 lsol = solvesdp(fcstr,z, get_yalmip_linear_options());
611 iF1 =
double(pFexp / pexp);
613 elseif lsol.problem == 1 || lsol.problem == 12
614 fprintf('Fitting MAMAP(2,2) F+S: linear program
is infeasible\n');
618 save(fname,'
map','p','F','S');
619 error('Fitting MAMAP(2,2) F+S: linear program error: %s, input saved to %s\n', lsol.info, fname);
625 % if feasible, solve nonlinear program
627 % declare variables for optimization in characteristic space
630 % initialize variables
631 if use_linprog_prepass
635 % define objective function
636 cobj = make_objective(vF1, vS11);
640 [cstr,~] = make_constraints_can1(p(1), vF1, vS11);
642 [cstr,~] = make_constraints_can2(p(1), vF1, vS11);
646 [~,cstr] = make_constraints_can1(p(1), vF1, vS11);
648 [~,cstr] = make_constraints_can2(p(1), vF1, vS11);
651 % solve nonlinear problem
652 csol = solvesdp(cstr, cobj, yalmip_nonlinear_opt);
656 save(fname,'
map','p','F','S');
657 error('Fitting MAMAP(2,2) F+S: solver (characteristic space) error: %s, input saved to %s\n', csol.info, fname);
662 fprintf('Fitting MAMAP(2,2) F+S: solver (characteristic space) objective = %e\n',
double(cobj));
665 % infeasible subproblem
672 function obj = make_objective(vF1, vS11)
673 obj = classWeights(1) * fsWeights(1) * (vF1/F(1) - 1)^2 + ...
674 classWeights(1) * fsWeights(2) * (vS11/S(1,1) - 1)^2;
677 function [q1n,q1d,q2n,q2d,q3] = make_q_can1(p, vF1, vS11)
678 % define coefficients of polynomials
679 q1n_F = p^2*(r1 - 1)*(r1*r2 - r2 + 1);
680 q1n_s = -r1*(h1 - h2 + h2*r1);
681 q1n_1 = h1*p^2*(r1*r2 - r2 + 1);
682 q1d_F = p*(r1 - 1)*(r1*r2 - r2 + 1);
683 q1d_1 = -p*(r1 - 1)*(h1 - h1*r2 + h2*r1);
684 q2n_Fsq = -p^2*(r1*r2 - r2 + 1)^2;
685 q2n_F = p^2*(r1*r2 - r2 + 1)*(2*h1 - h1*r1 - 2*h1*r2 + 3*h2*r1 - h2*r1^2 + h2*r1^2*r2 + h1*r1*r2 - h2*r1*r2);
686 q2n_s = -r1*(r2 - 1)*(h1 - h2 + h2*r1)^2;
687 q2n_1 = -p^2*(r1*r2 - r2 + 1)*(h2^2*r1 - h1^2*r2 + h1^2 + h1*h2*r1 - h1*h2*r1*r2);
688 q2d_F = p*r1*(r2 - 1)*(r1*r2 - r2 + 1)*(h1 - h2 + h2*r1);
689 q2d_1 = -p*r1*(r2 - 1)*(h1 - h1*r2 + h2*r1)*(h1 - h2 + h2*r1);
690 q3_F = -(p*(r1*r2 - r2 + 1))/(r1*r2*(h1 + h2*(r1 - 1)));
691 q3_1 = (p*(h1 + h2*r1)*(r1*r2 - r2 + 1))/(r1*r2*(h1 - h2 + h2*r1));
693 q1n = q1n_F * vF1 + q1n_s * vS11 + q1n_1;
694 q1d = q1d_F * vF1 + q1d_1;
695 q2n = q2n_Fsq * vF1^2 + q2n_F * vF1 + q2n_s * vS11 + q2n_1;
696 q2d = q2d_F * vF1 + q2d_1;
697 q3 = q3_F * vF1 + q3_1;
700 function [lcstr,rcstr] = make_constraints_can1(p, vF1, vS11)
702 % get fitting expressions
703 [q1n,q1d,q2n,q2d,q3] = make_q_can1(p, vF1, vS11);
705 % these are in common for the two subproblems
710 %scale = 1/abs(q2d{M1-tol});
712 % left constraints: F(1) < M1
716 %
do not change direction of inequalities
for q2
717 q2lb = scale*(q2n) >= 0;
718 q2ub = scale*(q2n) <= scale*(q2d);
720 % change direction of inequalities
for q2
721 q2lb = scale*(q2n) <= 0;
722 q2ub = scale*(q2n) >= scale*(q2d);
724 F1bound = vF1 <= M1-tol;
725 lcstr = [q1lb,q1ub,q2lb,q2ub,q3lb,q3ub,F1bound,vS11 >= 0, vS11 <= 1];
727 % right constraints: F(1) > M1
731 %
do not change direction of inequalities
for q2
732 q2lb = scale*(q2n) >= 0;
733 q2ub = scale*(q2n) <= scale*(q2d);
735 % change direction of inequalities
for q2
736 q2lb = scale*(q2n) <= 0;
737 q2ub = scale*(q2n) >= scale*(q2d);
739 F1bound = vF1 >= M1+tol;
740 rcstr = [q1lb,q1ub,q2lb,q2ub,q3lb,q3ub,F1bound,vS11 >= 0, vS11 <= 1];
743 function [q1n,q1d,q2,q3n,q3d] = make_q_can2(p, vF1, vS11)
744 % define coefficients of polynomials
745 q1n_F = p^2*(r1 - 1)*(r1 + r2 - r1*r2 - 2);
746 q1n_s = h1 - h2 + h2*r1;
747 q1n_1 = h1*p^2*(r1 + r2 - r1*r2 - 2);
748 q1d_F = p*(r1 - 1)*(r1 + r2 - r1*r2 - 2);
749 q1d_1 = p*(r1 - 1)*(h1 + h2 - h1*r2);
750 q2_F = p*(r1 + r2 - r1*r2 - 2)/(h2 - h1 + h1*r2 - h2*r1 - h2*r2 + h2*r1*r2);
751 q2_1 = -h2*p*(r1 + r2 - r1*r2 - 2)/(h2 - h1 + h1*r2 - h2*r1 - h2*r2 + h2*r1*r2);
752 q3n_Fsq = p^2*(r1 + r2 - r1*r2 - 2)^2;
753 q3n_F = p^2*(r1 + r2 - r1*r2 - 2)*(2*h1 + 2*h2 - h1*r2 - h2*r2 + h2*r1*r2);
754 q3n_s = -(r2 - 1)*(h1 - h2 + h2*r1)^2;
755 q3n_1 = -h2*p^2*(2*h1 - h1*r2 + h2*r1)*(r1 + r2 - r1*r2 - 2);
756 q3d_F = p*r2*(h1 - h2 + h2*r1)*(r1 + r2 - r1*r2 - 2);
757 q3d_1 = p*r2*(h1 + h2 - h1*r2)*(h1 - h2 + h2*r1);
759 q1n = q1n_F * vF1 + q1n_s * vS11 + q1n_1;
760 q1d = q1d_F * vF1 + q1d_1;
761 q2 = q2_F * vF1 + q2_1;
762 q3n = q3n_Fsq * vF1^2 + q3n_F * vF1 + q3n_s * vS11 + q3n_1;
763 q3d = q3d_F * vF1 + q3d_1;
766 function [lcstr,rcstr] = make_constraints_can2(p, vF1, vS11)
768 % get fitting expressions
769 [q1n,q1d,q2,q3n,q3d] = make_q_can2(p, vF1, vS11);
771 % these are in common
for the two subproblems
776 %scale = 1/abs(q3d{M1-tol});
778 % left constraints: F1 < M1
781 if (h1 - h2 + h2*r1) < 0
782 q3lb = scale*(q3n) <= 0;
783 q3ub = scale*(q3n) >= scale*(q3d);
785 q3lb = scale*(q3n) >= 0;
786 q3ub = scale*(q3n) <= scale*(q3d);
788 F1bound = vF1 <= M1-tol;
789 lcstr = [q1lb,q1ub,q2lb,q2ub,q3lb,q3ub,F1bound,vS11 >= 0, vS11 <= 1];
791 % rigt costraints: F1 > M1
794 if (h1 - h2 + h2*r1) < 0
795 q3lb = scale*(q3n) >= 0;
796 q3ub = scale*(q3n) <= scale*(q3d);
798 q3lb = scale*(q3n) <= 0;
799 q3ub = scale*(q3n) >= scale*(q3d);
801 F1bound = vF1 >= M1+tol;
802 rcstr = [q1lb,q1ub,q2lb,q2ub,q3lb,q3ub,F1bound,vS11 >= 0, vS11 <= 1];
805% used
for degenerate
case
806 function x = solve_quadprog()
807 fprintf('Fitting MAMAP(2,2) F+S: running quadratic programming solver...\n');
808 options = optimset('Algorithm','interior-point-convex','Display','none');
809 %[x,fx,xflag] = quadprog(H, h, A, b, [], [], [], [], [], options);
810 lb = 1e-6*ones( size(A,2),1);
811 ub = 1e6*ones( size(A,2),1);
812 [x,fx,xflag]=QP(H, h, A, b, Aeq, beq, lb, ub, options);
815 error('Quadratic programming solver failed: %d\n', exit);
817 fit_error = fx + length(x);
818 fprintf('Fitting MAMAP(2,2) F+S: error = %f\n', fit_error);