LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
qrf_rsrd.m
1function [result, x, fval, exitflag] = qrf_rsrd(params, objective, sense)
2% QRF_RSRD - Quadratic Reduction Framework for RS-RD blocking networks
3%
4% MATLAB port of the AMPL model qrboundsrsrd_skel.mod
5%
6% Usage:
7% [result, x, fval, exitflag] = qrf_rsrd(params)
8% [result, x, fval, exitflag] = qrf_rsrd(params, objective)
9% [result, x, fval, exitflag] = qrf_rsrd(params, objective, sense)
10%
11% Inputs:
12% params - Structure with model parameters:
13% .M - Number of queues
14% .N - Total population
15% .F - [M x 1] Capacity of each queue
16% .K - [M x 1] Number of phases for each queue
17% .mu - {M x 1} cell, each mu{i} is K(i) x K(i) completion rates
18% .v - {M x 1} cell, each v{i} is K(i) x K(i) background rates
19% .r - [M x M] Routing probabilities
20% .alpha - (optional) {M x 1} cell, each alpha{i} is [N x 1] load-dependent rates
21%
22% objective - (optional) 'U1min' (default), 'U1max', or queue index 1..M
23% sense - (optional) 'min' (default) or 'max'
24%
25% Outputs:
26% result - Structure with results:
27% .U - [M x 1] Utilization of each queue
28% .Ueff - [M x 1] Effective utilization
29% .pb - [M x 1] Blocking probability
30% .p2 - Decision variable tensor (marginal probabilities)
31% x - Raw solution vector
32% fval - Objective function value
33% exitflag - Solver exit flag
34
35 if nargin < 2 || isempty(objective)
36 objective = 'U1min';
37 end
38 if nargin < 3 || isempty(sense)
39 sense = 'min';
40 end
41
42 % Extract parameters
43 M = params.M;
44 N = params.N;
45 F = params.F(:);
46 K = params.K(:);
47 mu = params.mu;
48 v = params.v;
49 r = params.r;
50
51 % Default alpha (load-independent)
52 if isfield(params, 'alpha') && ~isempty(params.alpha)
53 alpha = params.alpha;
54 else
55 alpha = cell(M, 1);
56 for i = 1:M
57 alpha{i} = ones(N, 1);
58 end
59 end
60
61 % Compute transition rates q(i,j,k,h,n)
62 % q{i,j} is a K(i) x K(i) x (N+1) array
63 q = cell(M, M);
64 for i = 1:M
65 for j = 1:M
66 q{i,j} = zeros(K(i), K(i), N+1);
67 for ki = 1:K(i)
68 for hi = 1:K(i)
69 for n = 1:N % n=0 gives q=0
70 if j ~= i
71 q{i,j}(ki, hi, n+1) = r(i,j) * mu{i}(ki, hi) * alpha{i}(n);
72 else
73 q{i,j}(ki, hi, n+1) = (v{i}(ki, hi) + r(i,i) * mu{i}(ki, hi)) * alpha{i}(n);
74 end
75 end
76 end
77 end
78 end
79 end
80
81 %% Build variable indexing
82 % p2(j, nj, kj, i, ni, hi) for j in 1:M, nj in 0:F(j), kj in 1:K(j),
83 % i in 1:M, ni in 0:F(i), hi in 1:K(i)
84
85 % Count variables and build index map
86 fprintf('Building variable index map...\n');
87 varCount = 0;
88 p2idx = cell(M, 1);
89 for j = 1:M
90 p2idx{j} = cell(F(j)+1, K(j), M);
91 for nj = 0:F(j)
92 for kj = 1:K(j)
93 for i = 1:M
94 p2idx{j}{nj+1, kj, i} = zeros(F(i)+1, K(i));
95 for ni = 0:F(i)
96 for hi = 1:K(i)
97 varCount = varCount + 1;
98 p2idx{j}{nj+1, kj, i}(ni+1, hi) = varCount;
99 end
100 end
101 end
102 end
103 end
104 end
105 nP2Vars = varCount;
106
107 % U and Ueff variables: U(i,k,ni) and Ueff(i,k,ni) for ni >= 1
108 Uidx = cell(M, 1);
109 Ueffidx = cell(M, 1);
110 for i = 1:M
111 Uidx{i} = zeros(K(i), F(i));
112 Ueffidx{i} = zeros(K(i), F(i));
113 for ki = 1:K(i)
114 for ni = 1:F(i)
115 varCount = varCount + 1;
116 Uidx{i}(ki, ni) = varCount;
117 varCount = varCount + 1;
118 Ueffidx{i}(ki, ni) = varCount;
119 end
120 end
121 end
122
123 nVars = varCount;
124 fprintf('Total variables: %d (p2: %d, U/Ueff: %d)\n', nVars, nP2Vars, nVars - nP2Vars);
125
126 %% Build constraints
127 Aeq = [];
128 beq = [];
129 Aineq = [];
130 bineq = [];
131
132 % Helper to get variable index
133 getIdx = @(j, nj, kj, i, ni, hi) p2idx{j}{nj+1, kj, i}(ni+1, hi);
134 getUIdx = @(i, ki, ni) Uidx{i}(ki, ni);
135 getUeffIdx = @(i, ki, ni) Ueffidx{i}(ki, ni);
136
137 fprintf('Building constraints...\n');
138
139 %% ONE: Normalization - sum over nj, kj equals 1 for each j
140 fprintf(' ONE constraints...\n');
141 for j = 1:M
142 row = zeros(1, nVars);
143 for nj = 0:F(j)
144 for kj = 1:K(j)
145 idx = getIdx(j, nj, kj, j, nj, kj);
146 row(idx) = 1;
147 end
148 end
149 Aeq = [Aeq; row];
150 beq = [beq; 1];
151 end
152
153 %% ZERO constraints - fix infeasible states to zero
154 fprintf(' ZERO constraints...\n');
155 lb = zeros(nVars, 1);
156 ub = ones(nVars, 1);
157 % U and Ueff bounds are [0, 1]
158 for i = 1:M
159 for ki = 1:K(i)
160 for ni = 1:F(i)
161 ub(getUIdx(i, ki, ni)) = 1;
162 ub(getUeffIdx(i, ki, ni)) = 1;
163 end
164 end
165 end
166
167 for j = 1:M
168 for nj = 0:F(j)
169 for kj = 1:K(j)
170 for i = 1:M
171 for ni = 0:F(i)
172 for hi = 1:K(i)
173 idx = getIdx(j, nj, kj, i, ni, hi);
174
175 % ZERO1: i==j, nj==ni, h<>k
176 if i == j && nj == ni && hi ~= kj
177 ub(idx) = 0;
178 end
179
180 % ZERO2: i==j, nj<>ni
181 if i == j && nj ~= ni
182 ub(idx) = 0;
183 end
184
185 % ZERO3: i<>j, nj+ni > N
186 if i ~= j && nj + ni > N
187 ub(idx) = 0;
188 end
189
190 % ZERO6: i<>j, N-nj-ni > sum of other capacities
191 if i ~= j
192 sumOtherF = sum(F) - F(i) - F(j);
193 if N - nj - ni > sumOtherF
194 ub(idx) = 0;
195 end
196 end
197 end
198 end
199 end
200 end
201 end
202
203 % ZERO7: N-nj > sum of other capacities (for diagonal)
204 for nj = 0:F(j)
205 for kj = 1:K(j)
206 sumOtherF = sum(F) - F(j);
207 if N - nj > sumOtherF
208 idx = getIdx(j, nj, kj, j, nj, kj);
209 ub(idx) = 0;
210 end
211 end
212 end
213 end
214
215 %% SYMMETRY: p2(i,ni,hi,j,nj,kj) = p2(j,nj,kj,i,ni,hi)
216 fprintf(' SYMMETRY constraints...\n');
217 for j = 1:M
218 for nj = 0:F(j)
219 for kj = 1:K(j)
220 for i = 1:M
221 if i <= j
222 continue; % avoid duplicate constraints
223 end
224 for ni = 0:F(i)
225 for hi = 1:K(i)
226 idx1 = getIdx(j, nj, kj, i, ni, hi);
227 idx2 = getIdx(i, ni, hi, j, nj, kj);
228 if idx1 ~= idx2
229 row = zeros(1, nVars);
230 row(idx1) = 1;
231 row(idx2) = -1;
232 Aeq = [Aeq; row];
233 beq = [beq; 0];
234 end
235 end
236 end
237 end
238 end
239 end
240 end
241
242 %% MARGINALS: p2(j,nj,kj,j,nj,kj) = sum over ni,hi of p2(j,nj,kj,i,ni,hi)
243 fprintf(' MARGINALS constraints...\n');
244 for j = 1:M
245 for kj = 1:K(j)
246 for nj = 0:F(j)
247 for i = 1:M
248 if i == j
249 continue;
250 end
251 row = zeros(1, nVars);
252 idx_diag = getIdx(j, nj, kj, j, nj, kj);
253 row(idx_diag) = 1;
254 for ni = 0:F(i)
255 for hi = 1:K(i)
256 idx = getIdx(j, nj, kj, i, ni, hi);
257 row(idx) = row(idx) - 1;
258 end
259 end
260 Aeq = [Aeq; row];
261 beq = [beq; 0];
262 end
263 end
264 end
265 end
266
267 %% UCLASSIC: U(i,k,ni) = p2(i,ni,k,i,ni,k)
268 fprintf(' UCLASSIC constraints...\n');
269 for i = 1:M
270 for ki = 1:K(i)
271 for ni = 1:F(i)
272 row = zeros(1, nVars);
273 idx_U = getUIdx(i, ki, ni);
274 idx_p2 = getIdx(i, ni, ki, i, ni, ki);
275 row(idx_U) = 1;
276 row(idx_p2) = -1;
277 Aeq = [Aeq; row];
278 beq = [beq; 0];
279 end
280 end
281 end
282
283 %% UEFFS: Ueff(i,k,ni) = p2(i,ni,k,i,ni,k) - sum_{j: r(i,j)>0} r(i,j)*p2(i,ni,k,j,F(j),h)
284 fprintf(' UEFFS constraints...\n');
285 for i = 1:M
286 for ki = 1:K(i)
287 for ni = 1:F(i)
288 row = zeros(1, nVars);
289 idx_Ueff = getUeffIdx(i, ki, ni);
290 idx_p2_diag = getIdx(i, ni, ki, i, ni, ki);
291 row(idx_Ueff) = 1;
292 row(idx_p2_diag) = -1;
293 % Add blocking terms
294 for j = 1:M
295 if j ~= i && r(i,j) > 0
296 for hj = 1:K(j)
297 idx_block = getIdx(i, ni, ki, j, F(j), hj);
298 row(idx_block) = r(i,j);
299 end
300 end
301 end
302 Aeq = [Aeq; row];
303 beq = [beq; 0];
304 end
305 end
306 end
307
308 %% THM2: Phase balance (Theorem 1 - THM:sdeffective)
309 % sum_{ni} (sum_{j<>i, h<>k} q(i,j,k,h,ni)*Ueff(i,k,ni) + sum_{h<>k} q(i,i,k,h,ni)*U(i,k,ni))
310 % = sum_{ni} (sum_{j<>i, h<>k} q(i,j,h,k,ni)*Ueff(i,h,ni) + sum_{h<>k} q(i,i,h,k,ni)*U(i,h,ni))
311 fprintf(' THM2 (Phase balance) constraints...\n');
312 for i = 1:M
313 for ki = 1:K(i)
314 row = zeros(1, nVars);
315 % LHS terms
316 for ni = 1:F(i)
317 % Terms with Ueff (j<>i)
318 for j = 1:M
319 if j == i
320 continue;
321 end
322 for hi = 1:K(i)
323 if hi == ki
324 continue;
325 end
326 coef = q{i,j}(ki, hi, ni+1);
327 idx_Ueff = getUeffIdx(i, ki, ni);
328 row(idx_Ueff) = row(idx_Ueff) + coef;
329 end
330 end
331 % Terms with U (j==i, h<>k)
332 for hi = 1:K(i)
333 if hi == ki
334 continue;
335 end
336 coef = q{i,i}(ki, hi, ni+1);
337 idx_p2 = getIdx(i, ni, ki, i, ni, ki);
338 row(idx_p2) = row(idx_p2) + coef;
339 end
340 end
341 % RHS terms (subtract)
342 for ni = 1:F(i)
343 % Terms with Ueff (j<>i)
344 for j = 1:M
345 if j == i
346 continue;
347 end
348 for hi = 1:K(i)
349 if hi == ki
350 continue;
351 end
352 coef = q{i,j}(hi, ki, ni+1);
353 idx_Ueff = getUeffIdx(i, hi, ni);
354 row(idx_Ueff) = row(idx_Ueff) - coef;
355 end
356 end
357 % Terms with U (j==i, h<>k)
358 for hi = 1:K(i)
359 if hi == ki
360 continue;
361 end
362 coef = q{i,i}(hi, ki, ni+1);
363 idx_p2 = getIdx(i, ni, hi, i, ni, hi);
364 row(idx_p2) = row(idx_p2) - coef;
365 end
366 end
367 Aeq = [Aeq; row];
368 beq = [beq; 0];
369 end
370 end
371
372 %% THM1: Population constraint (Theorem 2)
373 % sum_i sum_ni sum_hi ni * p2(j,nj,kj,i,ni,hi) = N * p2(j,nj,kj,j,nj,kj)
374 fprintf(' THM1 (Population) constraints...\n');
375 for j = 1:M
376 for kj = 1:K(j)
377 for nj = 0:F(j)
378 row = zeros(1, nVars);
379 % RHS: -N * p2(j,nj,kj,j,nj,kj)
380 idx_diag = getIdx(j, nj, kj, j, nj, kj);
381 row(idx_diag) = -N;
382 % LHS: sum
383 for i = 1:M
384 for ni = 1:F(i) % ni >= 1
385 for hi = 1:K(i)
386 idx = getIdx(j, nj, kj, i, ni, hi);
387 row(idx) = row(idx) + ni;
388 end
389 end
390 end
391 Aeq = [Aeq; row];
392 beq = [beq; 0];
393 end
394 end
395 end
396
397 %% THM3a: Marginal balance for ni in 1:F(i)-1
398 fprintf(' THM3a (Marginal balance) constraints...\n');
399 for i = 1:M
400 for ni = 1:(F(i)-1)
401 row = zeros(1, nVars);
402 % LHS: arrivals to queue i
403 for j = 1:M
404 if j == i
405 continue;
406 end
407 for kj = 1:K(j)
408 for hj = 1:K(j)
409 for ui = 1:K(i)
410 for nj = 1:F(j)
411 idx = getIdx(j, nj, kj, i, ni, ui);
412 row(idx) = row(idx) + q{j,i}(kj, hj, nj+1);
413 end
414 end
415 end
416 end
417 end
418 % RHS: departures from queue i at ni+1
419 for j = 1:M
420 if j == i
421 continue;
422 end
423 for ki = 1:K(i)
424 for hi = 1:K(i)
425 for uj = 1:K(j)
426 for nj = 0:(F(j)-1)
427 idx = getIdx(i, ni+1, ki, j, nj, uj);
428 row(idx) = row(idx) - q{i,j}(ki, hi, ni+2);
429 end
430 end
431 end
432 end
433 end
434 Aeq = [Aeq; row];
435 beq = [beq; 0];
436 end
437 end
438
439 %% THM3b: Marginal balance for ni=0, per phase
440 fprintf(' THM3b (Marginal balance ni=0) constraints...\n');
441 for i = 1:M
442 for ui = 1:K(i)
443 row = zeros(1, nVars);
444 % LHS: arrivals to queue i at ni=0
445 for j = 1:M
446 if j == i
447 continue;
448 end
449 for kj = 1:K(j)
450 for hj = 1:K(j)
451 for nj = 1:F(j)
452 idx = getIdx(j, nj, kj, i, 0, ui);
453 row(idx) = row(idx) + q{j,i}(kj, hj, nj+1);
454 end
455 end
456 end
457 end
458 % RHS: departures from queue i at ni=1
459 for j = 1:M
460 if j == i
461 continue;
462 end
463 for ki = 1:K(i)
464 for nj = 0:(F(j)-1)
465 for hj = 1:K(j)
466 idx = getIdx(i, 1, ki, j, nj, hj);
467 row(idx) = row(idx) - q{i,j}(ki, ui, 2);
468 end
469 end
470 end
471 end
472 Aeq = [Aeq; row];
473 beq = [beq; 0];
474 end
475 end
476
477 %% QBAL: Queue balance constraint
478 % This is a complex balance equation that tightens the bounds
479 fprintf(' QBAL (Queue balance) constraints...\n');
480 for i = 1:M
481 for ki = 1:K(i)
482 row = zeros(1, nVars);
483
484 % LHS Term 1: sum{h<>k} sum{j<>i} sum{ni} sum{u} sum{nj} q[i,j,k,h,ni]*ni*p2[i,ni,k,j,nj,u]
485 for hi = 1:K(i)
486 if hi == ki
487 continue;
488 end
489 for j = 1:M
490 if j == i
491 continue;
492 end
493 for ni = 1:F(i)
494 for uj = 1:K(j)
495 for nj = 0:(F(j)-1)
496 coef = q{i,j}(ki, hi, ni+1) * ni;
497 idx = getIdx(i, ni, ki, j, nj, uj);
498 row(idx) = row(idx) + coef;
499 end
500 end
501 end
502 end
503 end
504
505 % LHS Term 2: sum{h<>k} sum{ni} q[i,i,k,h,ni]*ni*p2[i,ni,k,i,ni,k]
506 for hi = 1:K(i)
507 if hi == ki
508 continue;
509 end
510 for ni = 1:F(i)
511 coef = q{i,i}(ki, hi, ni+1) * ni;
512 idx = getIdx(i, ni, ki, i, ni, ki);
513 row(idx) = row(idx) + coef;
514 end
515 end
516
517 % LHS Term 3: sum{j<>i} sum{h} sum{ni} sum{u} sum{nj<=min(F(j)-1,N-ni)} q[i,j,h,k,ni]*p2[i,ni,h,j,nj,u]
518 for j = 1:M
519 if j == i
520 continue;
521 end
522 for hi = 1:K(i)
523 for ni = 1:F(i)
524 for uj = 1:K(j)
525 maxNj = min(F(j)-1, N-ni);
526 for nj = 0:maxNj
527 coef = q{i,j}(hi, ki, ni+1);
528 idx = getIdx(i, ni, hi, j, nj, uj);
529 row(idx) = row(idx) + coef;
530 end
531 end
532 end
533 end
534 end
535
536 % RHS Term 1: sum{j<>i} sum{h} sum{ni<=F(i)-1} sum{u} sum{nj>=1} q[j,i,h,u,nj]*p2[i,ni,k,j,nj,h]
537 for j = 1:M
538 if j == i
539 continue;
540 end
541 for hj = 1:K(j)
542 for ni = 0:(F(i)-1)
543 for uj = 1:K(j)
544 for nj = 1:F(j)
545 coef = q{j,i}(hj, uj, nj+1);
546 idx = getIdx(i, ni, ki, j, nj, hj);
547 row(idx) = row(idx) - coef;
548 end
549 end
550 end
551 end
552 end
553
554 % RHS Term 2: sum{h<>k} sum{ni} q[i,i,h,k,ni]*ni*p2[i,ni,h,i,ni,h]
555 for hi = 1:K(i)
556 if hi == ki
557 continue;
558 end
559 for ni = 1:F(i)
560 coef = q{i,i}(hi, ki, ni+1) * ni;
561 idx = getIdx(i, ni, hi, i, ni, hi);
562 row(idx) = row(idx) - coef;
563 end
564 end
565
566 % RHS Term 3: sum{h<>k} sum{j<>i} sum{ni} sum{u} sum{nj} q[i,j,h,k,ni]*ni*p2[i,ni,h,j,nj,u]
567 for hi = 1:K(i)
568 if hi == ki
569 continue;
570 end
571 for j = 1:M
572 if j == i
573 continue;
574 end
575 for ni = 1:F(i)
576 for uj = 1:K(j)
577 for nj = 0:(F(j)-1)
578 coef = q{i,j}(hi, ki, ni+1) * ni;
579 idx = getIdx(i, ni, hi, j, nj, uj);
580 row(idx) = row(idx) - coef;
581 end
582 end
583 end
584 end
585 end
586
587 Aeq = [Aeq; row];
588 beq = [beq; 0];
589 end
590 end
591
592 %% THM4: Queue-length bound inequality (Theorem 5)
593 % sum_t sum_ht sum_nj sum_nt nt*p2(j,nj,kj,t,nt,ht) >= N * sum_hi sum_nj sum_ni p2(j,nj,kj,i,ni,hi)
594 fprintf(' THM4 (Queue-length bound) constraints...\n');
595 for j = 1:M
596 for kj = 1:K(j)
597 for i = 1:M
598 row = zeros(1, nVars);
599 % LHS: sum_t sum_ht sum_nj sum_nt nt * p2
600 for t = 1:M
601 for ht = 1:K(t)
602 for nj = 0:F(j)
603 for nt = 1:F(t) % nt >= 1
604 idx = getIdx(j, nj, kj, t, nt, ht);
605 row(idx) = row(idx) + nt;
606 end
607 end
608 end
609 end
610 % RHS: -N * sum
611 for hi = 1:K(i)
612 for nj = 0:F(j)
613 for ni = 1:F(i) % ni >= 1
614 idx = getIdx(j, nj, kj, i, ni, hi);
615 row(idx) = row(idx) - N;
616 end
617 end
618 end
619 Aineq = [Aineq; -row]; % >= becomes <= with negation
620 bineq = [bineq; 0];
621 end
622 end
623 end
624
625 %% Build objective function
626 fprintf('Building objective function...\n');
627 c = zeros(nVars, 1);
628
629 if ischar(objective)
630 if strcmp(objective, 'U1min') || strcmp(objective, 'U1max')
631 targetQueue = 1;
632 else
633 error('Unknown objective: %s', objective);
634 end
635 else
636 targetQueue = objective;
637 end
638
639 % Utilization = sum over k, n of p2(i,n,k,i,n,k)
640 for ki = 1:K(targetQueue)
641 for ni = 1:F(targetQueue)
642 idx = getIdx(targetQueue, ni, ki, targetQueue, ni, ki);
643 c(idx) = 1;
644 end
645 end
646
647 if strcmp(sense, 'max')
648 c = -c;
649 end
650
651 %% Solve LP
652 fprintf('Solving LP with %d variables and %d equality + %d inequality constraints...\n', ...
653 nVars, size(Aeq, 1), size(Aineq, 1));
654
655 options = optimoptions('linprog', 'Display', 'final', 'Algorithm', 'interior-point');
656
657 [x, fval, exitflag] = linprog(c, Aineq, bineq, Aeq, beq, lb, ub, options);
658
659 if strcmp(sense, 'max')
660 fval = -fval;
661 end
662
663 %% Extract results
664 result = struct();
665 result.objective = fval;
666 result.exitflag = exitflag;
667
668 if exitflag > 0
669 % Compute utilizations from U and Ueff variables
670 result.U = zeros(M, 1);
671 result.Ueff = zeros(M, 1);
672 result.pb = zeros(M, 1);
673
674 for i = 1:M
675 % U(i) = sum over k, ni of U(i,k,ni)
676 for ki = 1:K(i)
677 for ni = 1:F(i)
678 idx_U = getUIdx(i, ki, ni);
679 idx_Ueff = getUeffIdx(i, ki, ni);
680 result.U(i) = result.U(i) + x(idx_U);
681 result.Ueff(i) = result.Ueff(i) + x(idx_Ueff);
682 end
683 end
684 result.pb(i) = result.U(i) - result.Ueff(i);
685 end
686
687 % Store p2 as a function handle for easy access
688 result.getP2 = @(j, nj, kj, i, ni, hi) x(getIdx(j, nj, kj, i, ni, hi));
689 end
690
691 fprintf('\n=== Results ===\n');
692 fprintf('Objective value: %f\n', fval);
693 fprintf('Exit flag: %d\n', exitflag);
694 if exitflag > 0
695 fprintf('\nUtilizations:\n');
696 for i = 1:M
697 fprintf(' Queue %d: U = %.6f, Ueff = %.6f, pb = %.6f\n', ...
698 i, result.U(i), result.Ueff(i), result.pb(i));
699 end
700 end
701end