LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
mapqn_bnd_qr_ld.m
1function [result, x, fval, exitflag] = mapqn_bnd_qr_ld(params, objective_queue, objective_phase, objective_n, sense)
2% MAPQN_BND_QR_LD - Quadratic Reduction Bounds for Load-Dependent Systems
3%
4% MATLAB port of bnd_qr_ld.py (Quadratic Reduction Bounds for Load-Dependent
5% MAP queueing networks with load-dependent service rates).
6%
7% Usage:
8% [result, x, fval, exitflag] = mapqn_bnd_qr_ld(params, objective_queue, objective_phase, objective_n)
9% [result, x, fval, exitflag] = mapqn_bnd_qr_ld(params, objective_queue, objective_phase, objective_n, sense)
10%
11% Inputs:
12% params - Structure with model parameters:
13% .M - Number of queues
14% .N - Total population
15% .K - [M x 1] Number of phases for each queue
16% .mu - {M x 1} cell, each mu{i} is K(i) x K(i) completion rates
17% .v - {M x 1} cell, each v{i} is K(i) x K(i) background rates
18% .alpha - [M x N] load-dependent rates
19% .r - [M x M] Routing probabilities
20% .verbose - (optional) boolean, default true
21%
22% objective_queue - Queue index for objective (1-based, 1..M)
23% objective_phase - Phase index for objective (1-based, 1..K(objective_queue))
24% objective_n - Population level for objective (0..N)
25% sense - (optional) 'max' (default) or 'min'
26%
27% Outputs:
28% result - Structure with results:
29% .objective - Objective function value
30% .exitflag - Solver exit flag
31% .p2marginals - Marginal probabilities
32% x - Raw solution vector
33% fval - Objective function value
34% exitflag - Solver exit flag
35
36 if nargin < 5 || isempty(sense)
37 sense = 'max';
38 end
39
40 % Extract parameters
41 M = params.M;
42 N = params.N;
43 K = params.K(:);
44 mu = params.mu;
45 v = params.v;
46 alpha = params.alpha;
47 r = params.r;
48 if isfield(params, 'verbose')
49 verbose = params.verbose;
50 else
51 verbose = true;
52 end
53
54 % Validate objective indices
55 if objective_queue < 1 || objective_queue > M
56 error('objective_queue must be in range 1..%d', M);
57 end
58 if objective_phase < 1 || objective_phase > K(objective_queue)
59 error('objective_phase must be in range 1..%d', K(objective_queue));
60 end
61 if objective_n < 0 || objective_n > N
62 error('objective_n must be in range 0..%d', N);
63 end
64
65 %% Define q function (1-based indices, load-dependent)
66 % q_val(i,j,k,h,n): transition rate from queue i phase k to queue j phase h
67 % when queue i has n customers
68 q_val = @(i, j, k, h, n) q_func(i, j, k, h, n, mu, v, alpha, r);
69
70 %% Build variable indexing
71 % p2(j, nj, k, i, ni, h): joint probabilities
72 % j in 1:M, nj in 0:N, k in 1:K(j), i in 1:M, ni in 0:N, h in 1:K(i)
73 if verbose; fprintf('Building variable index map...\n'); end
74 varCount = 0;
75
76 % p2 variables
77 p2idx = cell(M, 1);
78 for j = 1:M
79 p2idx{j} = cell(N+1, K(j), M);
80 for nj = 0:N
81 for kj = 1:K(j)
82 for i = 1:M
83 p2idx{j}{nj+1, kj, i} = zeros(N+1, K(i));
84 for ni = 0:N
85 for hi = 1:K(i)
86 varCount = varCount + 1;
87 p2idx{j}{nj+1, kj, i}(ni+1, hi) = varCount;
88 end
89 end
90 end
91 end
92 end
93 end
94
95 nVars = varCount;
96 if verbose; fprintf('Total variables: %d\n', nVars); end
97
98 %% Helper function
99 getP2Idx = @(j, nj, kj, i, ni, hi) p2idx{j}{nj+1, kj, i}(ni+1, hi);
100
101 %% Initialize bounds
102 lb = zeros(nVars, 1);
103 ub = ones(nVars, 1);
104
105 %% Build constraints
106 Aeq = [];
107 beq = [];
108 Aineq = [];
109 bineq = [];
110
111 if verbose; fprintf('Building constraints...\n'); end
112
113 %% ZERO constraints - fix infeasible states via ub=0
114 if verbose; fprintf(' ZERO constraints...\n'); end
115 for j = 1:M
116 for nj = 0:N
117 for kj = 1:K(j)
118 for i = 1:M
119 for ni = 0:N
120 for hi = 1:K(i)
121 idx = getP2Idx(j, nj, kj, i, ni, hi);
122
123 % ZERO1: i==j, nj==ni, h~=k
124 if i == j && nj == ni && hi ~= kj
125 ub(idx) = 0;
126 end
127
128 % ZERO2: i==j, nj~=ni
129 if i == j && nj ~= ni
130 ub(idx) = 0;
131 end
132
133 % ZERO3: i~=j, nj+ni > N
134 if i ~= j && nj + ni > N
135 ub(idx) = 0;
136 end
137 end
138 end
139 end
140 end
141 end
142 end
143
144 %% ONE: Normalization
145 % sum{nj,k} p2(j,nj,k,j,nj,k) = 1 for each j
146 if verbose; fprintf(' ONE constraints...\n'); end
147 for j = 1:M
148 row = zeros(1, nVars);
149 for nj = 0:N
150 for kj = 1:K(j)
151 idx = getP2Idx(j, nj, kj, j, nj, kj);
152 row(idx) = 1;
153 end
154 end
155 Aeq = [Aeq; row];
156 beq = [beq; 1];
157 end
158
159 %% SYMMETRY: p2(i,ni,h,j,nj,k) = p2(j,nj,k,i,ni,h) (only one ordering)
160 if verbose; fprintf(' SYMMETRY constraints...\n'); end
161 for j = 1:M
162 for nj = 0:N
163 for kj = 1:K(j)
164 for i = 1:M
165 if i <= j && ~(i == j)
166 % skip: we only do i > j for ordering
167 end
168 if i <= j
169 continue;
170 end
171 for ni = 0:N
172 if i ~= j && nj + ni > N
173 continue;
174 end
175 for hi = 1:K(i)
176 idx1 = getP2Idx(j, nj, kj, i, ni, hi);
177 idx2 = getP2Idx(i, ni, hi, j, nj, kj);
178 if ub(idx1) == 0 && ub(idx2) == 0
179 continue;
180 end
181 if idx1 ~= idx2
182 row = zeros(1, nVars);
183 row(idx1) = 1;
184 row(idx2) = -1;
185 Aeq = [Aeq; row];
186 beq = [beq; 0];
187 end
188 end
189 end
190 end
191 end
192 end
193 end
194
195 %% MARGINALS: p2(j,nj,k,j,nj,k) = sum{ni=0..N-nj, h} p2(j,nj,k,i,ni,h) for i~=j
196 if verbose; fprintf(' MARGINALS constraints...\n'); end
197 for j = 1:M
198 for kj = 1:K(j)
199 for nj = 0:N
200 for i = 1:M
201 if i == j
202 continue;
203 end
204 row = zeros(1, nVars);
205 idx_diag = getP2Idx(j, nj, kj, j, nj, kj);
206 row(idx_diag) = 1;
207 for ni = 0:(N - nj)
208 for hi = 1:K(i)
209 idx = getP2Idx(j, nj, kj, i, ni, hi);
210 row(idx) = row(idx) - 1;
211 end
212 end
213 Aeq = [Aeq; row];
214 beq = [beq; 0];
215 end
216 end
217 end
218 end
219
220 %% THM1: Queue length theorem
221 % for each (j,k): sum{i,nj>=1,ni>=1,h} ni*p2(j,nj,k,i,ni,h)
222 % - N*sum{nj>=1} p2(j,nj,k,j,nj,k) = 0
223 if verbose; fprintf(' THM1 (Queue length) constraints...\n'); end
224 for j = 1:M
225 for kj = 1:K(j)
226 row = zeros(1, nVars);
227 % LHS: sum over i, nj>=1, ni>=1, h of ni*p2
228 for i = 1:M
229 for nj = 1:N
230 for ni = 1:N
231 for hi = 1:K(i)
232 idx = getP2Idx(j, nj, kj, i, ni, hi);
233 row(idx) = row(idx) + ni;
234 end
235 end
236 end
237 end
238 % RHS: -N * sum{nj>=1} p2(j,nj,k,j,nj,k)
239 for nj = 1:N
240 idx = getP2Idx(j, nj, kj, j, nj, kj);
241 row(idx) = row(idx) - N;
242 end
243 Aeq = [Aeq; row];
244 beq = [beq; 0];
245 end
246 end
247
248 %% THM1c: Empty queue theorem
249 % for each (j,k): sum{i,ni>=1,h} ni*p2(j,0,k,i,ni,h) - N*p2(j,0,k,j,0,k) = 0
250 if verbose; fprintf(' THM1c (Empty queue) constraints...\n'); end
251 for j = 1:M
252 for kj = 1:K(j)
253 row = zeros(1, nVars);
254 % LHS
255 for i = 1:M
256 for ni = 1:N
257 for hi = 1:K(i)
258 idx = getP2Idx(j, 0, kj, i, ni, hi);
259 row(idx) = row(idx) + ni;
260 end
261 end
262 end
263 % RHS
264 idx = getP2Idx(j, 0, kj, j, 0, kj);
265 row(idx) = row(idx) - N;
266 Aeq = [Aeq; row];
267 beq = [beq; 0];
268 end
269 end
270
271 %% PC2: Second moment constraint
272 % sum{i,j,ni>=1,nj>=1,h,k} nj*ni*p2(j,nj,k,i,ni,h) = N^2
273 if verbose; fprintf(' PC2 (Second moment) constraint...\n'); end
274 row = zeros(1, nVars);
275 for i = 1:M
276 for j = 1:M
277 for ni = 1:N
278 for nj = 1:N
279 for hi = 1:K(i)
280 for kj = 1:K(j)
281 idx = getP2Idx(j, nj, kj, i, ni, hi);
282 row(idx) = row(idx) + nj * ni;
283 end
284 end
285 end
286 end
287 end
288 end
289 Aeq = [Aeq; row];
290 beq = [beq; N^2];
291
292 %% THM2: Phase balance
293 % for each (i,k): sum{j,h: j~=i or h~=k, ni>=1}
294 % (q(i,j,k,h,ni)*p2(i,ni,k,i,ni,k) - q(i,j,h,k,ni)*p2(i,ni,h,i,ni,h)) = 0
295 if verbose; fprintf(' THM2 (Phase balance) constraints...\n'); end
296 for i = 1:M
297 for ki = 1:K(i)
298 row = zeros(1, nVars);
299 for j = 1:M
300 for hi = 1:K(i)
301 if ~(hi == ki && j == i)
302 for ni = 1:N
303 q_out = q_val(i, j, ki, hi, ni);
304 q_in = q_val(i, j, hi, ki, ni);
305 idx_k = getP2Idx(i, ni, ki, i, ni, ki);
306 idx_h = getP2Idx(i, ni, hi, i, ni, hi);
307 row(idx_k) = row(idx_k) + q_out;
308 row(idx_h) = row(idx_h) - q_in;
309 end
310 end
311 end
312 end
313 Aeq = [Aeq; row];
314 beq = [beq; 0];
315 end
316 end
317
318 %% THM3a: Flow balance for ni = 1..N-1
319 % for each (i, ni in 1..N-1):
320 % sum{j~=i, k, h, u, nj=1..N-ni} q(j,i,k,h,nj)*p2(j,nj,k,i,ni,u)
321 % - sum{j~=i, k, h} q(i,j,k,h,ni+1)*p2(i,ni+1,k,i,ni+1,k) = 0
322 if verbose; fprintf(' THM3a (Flow balance ni=1..N-1) constraints...\n'); end
323 for i = 1:M
324 for ni = 1:(N-1)
325 row = zeros(1, nVars);
326 % Incoming flow
327 for j = 1:M
328 if j ~= i
329 for kj = 1:K(j)
330 for hj = 1:K(j)
331 for u = 1:K(i)
332 for nj = 1:(N - ni)
333 qv = q_val(j, i, kj, hj, nj);
334 idx = getP2Idx(j, nj, kj, i, ni, u);
335 row(idx) = row(idx) + qv;
336 end
337 end
338 end
339 end
340 end
341 end
342 % Outgoing flow
343 for j = 1:M
344 if j ~= i
345 for ki = 1:K(i)
346 for hi = 1:K(i)
347 qv = q_val(i, j, ki, hi, ni + 1);
348 idx = getP2Idx(i, ni + 1, ki, i, ni + 1, ki);
349 row(idx) = row(idx) - qv;
350 end
351 end
352 end
353 end
354 Aeq = [Aeq; row];
355 beq = [beq; 0];
356 end
357 end
358
359 %% THM3b: Flow balance for ni = 0
360 % for each (i, u):
361 % sum{j~=i, k, h, nj=1..N} q(j,i,k,h,nj)*p2(j,nj,k,i,0,u)
362 % - sum{j~=i, k} q(i,j,k,u,1)*p2(i,1,k,i,1,k) = 0
363 if verbose; fprintf(' THM3b (Flow balance ni=0) constraints...\n'); end
364 for i = 1:M
365 for u = 1:K(i)
366 row = zeros(1, nVars);
367 % Incoming to empty queue
368 for j = 1:M
369 if j ~= i
370 for kj = 1:K(j)
371 for hj = 1:K(j)
372 for nj = 1:N
373 qv = q_val(j, i, kj, hj, nj);
374 idx = getP2Idx(j, nj, kj, i, 0, u);
375 row(idx) = row(idx) + qv;
376 end
377 end
378 end
379 end
380 end
381 % Outgoing from single customer
382 for j = 1:M
383 if j ~= i
384 for ki = 1:K(i)
385 qv = q_val(i, j, ki, u, 1);
386 idx = getP2Idx(i, 1, ki, i, 1, ki);
387 row(idx) = row(idx) - qv;
388 end
389 end
390 end
391 Aeq = [Aeq; row];
392 beq = [beq; 0];
393 end
394 end
395
396 %% QBAL: Queue balance (from AMPL bnd_quadraticreduction_ld.mod)
397 % for each (i, k):
398 % LHS1: sum{h~=k, j, ni>=1} q(i,j,k,h,ni)*ni*p2(i,ni,k,i,ni,k)
399 % + LHS2: sum{j~=i, h, ni>=1} q(i,j,h,k,ni)*p2(i,ni,h,i,ni,h)
400 % = RHS1: sum{j~=i, u in K(j), w in K(j), nj>=1} q(j,i,u,w,nj)*(p2(j,nj,u,i,0,k) + sum{ni>=1} p2(i,ni,k,j,nj,u))
401 % + RHS2: sum{h~=k, j, ni>=1} q(i,j,h,k,ni)*ni*p2(i,ni,h,i,ni,h)
402 if verbose; fprintf(' QBAL (Queue balance) constraints...\n'); end
403 for i = 1:M
404 for ki = 1:K(i)
405 row = zeros(1, nVars);
406 % LHS1: sum{h~=k, j, ni>=1} q(i,j,k,h,ni)*ni*p2(i,ni,k,i,ni,k)
407 for hi = 1:K(i)
408 if hi ~= ki
409 for j = 1:M
410 for ni = 1:N
411 qv = q_val(i, j, ki, hi, ni);
412 idx = getP2Idx(i, ni, ki, i, ni, ki);
413 row(idx) = row(idx) + qv * ni;
414 end
415 end
416 end
417 end
418 % LHS2: sum{j~=i, h, ni>=1} q(i,j,h,k,ni)*p2(i,ni,h,i,ni,h)
419 for j = 1:M
420 if j ~= i
421 for hi = 1:K(i)
422 for ni = 1:N
423 qv = q_val(i, j, hi, ki, ni);
424 idx = getP2Idx(i, ni, hi, i, ni, hi);
425 row(idx) = row(idx) + qv;
426 end
427 end
428 end
429 end
430 % -RHS1 part a: sum{j~=i, u in K(j), w in K(j), nj>=1} q(j,i,u,w,nj)*p2(j,nj,u,i,0,k)
431 for j = 1:M
432 if j ~= i
433 for u = 1:K(j)
434 for w = 1:K(j)
435 for nj = 1:N
436 qv = q_val(j, i, u, w, nj);
437 idx = getP2Idx(j, nj, u, i, 0, ki);
438 row(idx) = row(idx) - qv;
439 end
440 end
441 end
442 end
443 end
444 % -RHS1 part b: sum{j~=i, u in K(j), w in K(j), nj>=1, ni>=1} q(j,i,u,w,nj)*p2(i,ni,k,j,nj,u)
445 for j = 1:M
446 if j ~= i
447 for u = 1:K(j)
448 for w = 1:K(j)
449 for nj = 1:N
450 qv = q_val(j, i, u, w, nj);
451 for ni = 1:N
452 idx = getP2Idx(i, ni, ki, j, nj, u);
453 row(idx) = row(idx) - qv;
454 end
455 end
456 end
457 end
458 end
459 end
460 % -RHS2: sum{h~=k, j, ni>=1} q(i,j,h,k,ni)*ni*p2(i,ni,h,i,ni,h)
461 for hi = 1:K(i)
462 if hi ~= ki
463 for j = 1:M
464 for ni = 1:N
465 qv = q_val(i, j, hi, ki, ni);
466 idx = getP2Idx(i, ni, hi, i, ni, hi);
467 row(idx) = row(idx) - qv * ni;
468 end
469 end
470 end
471 end
472 Aeq = [Aeq; row];
473 beq = [beq; 0];
474 end
475 end
476
477 %% COR1a: Correlation constraint (from AMPL bnd_quadraticreduction_ld.mod)
478 % for each (i, kstar, ni_c in 0..N-2)
479 if verbose; fprintf(' COR1a constraints...\n'); end
480 for i = 1:M
481 for kstar = 1:K(i)
482 for ni_c = 0:(N-2)
483 row = zeros(1, nVars);
484 % A: sum{j~=i,kj,hj,u~=kstar,nj=1..N-ni_c} q(j,i,kj,hj,nj)*p2(j,nj,kj,i,ni_c,u)
485 for j = 1:M
486 if j ~= i
487 for kj = 1:K(j)
488 for hj = 1:K(j)
489 for u = 1:K(i)
490 if u ~= kstar
491 for nj = 1:(N-ni_c)
492 qv = q_val(j, i, kj, hj, nj);
493 idx = getP2Idx(j, nj, kj, i, ni_c, u);
494 row(idx) = row(idx) + qv;
495 end
496 end
497 end
498 end
499 end
500 end
501 end
502 % B: sum{j~=i,kj,hj,nj=1..N-ni_c} q(j,i,kj,hj,nj)*p2(j,nj,kj,i,ni_c+1,kstar)
503 for j = 1:M
504 if j ~= i
505 for kj = 1:K(j)
506 for hj = 1:K(j)
507 for nj = 1:(N-ni_c)
508 qv = q_val(j, i, kj, hj, nj);
509 idx = getP2Idx(j, nj, kj, i, ni_c+1, kstar);
510 row(idx) = row(idx) + qv;
511 end
512 end
513 end
514 end
515 end
516 % C: sum{k~=kstar} q(i,i,kstar,k,ni_c+1)*p2(i,ni_c+1,kstar,i,ni_c+1,kstar)
517 for k2 = 1:K(i)
518 if k2 ~= kstar
519 qv = q_val(i, i, kstar, k2, ni_c+1);
520 idx = getP2Idx(i, ni_c+1, kstar, i, ni_c+1, kstar);
521 row(idx) = row(idx) + qv;
522 end
523 end
524 % -D: -sum{j~=i,k~=kstar} q(i,j,k,k,ni_c+1)*p2(i,ni_c+1,k,i,ni_c+1,k)
525 for j = 1:M
526 if j ~= i
527 for k2 = 1:K(i)
528 if k2 ~= kstar
529 qv = q_val(i, j, k2, k2, ni_c+1);
530 idx = getP2Idx(i, ni_c+1, k2, i, ni_c+1, k2);
531 row(idx) = row(idx) - qv;
532 end
533 end
534 end
535 end
536 % -E: -sum{j~=i,k~=kstar,h~=k} q(i,j,k,h,ni_c+1)*p2(i,ni_c+1,k,i,ni_c+1,k)
537 for j = 1:M
538 if j ~= i
539 for k2 = 1:K(i)
540 if k2 ~= kstar
541 for h2 = 1:K(i)
542 if h2 ~= k2
543 qv = q_val(i, j, k2, h2, ni_c+1);
544 idx = getP2Idx(i, ni_c+1, k2, i, ni_c+1, k2);
545 row(idx) = row(idx) - qv;
546 end
547 end
548 end
549 end
550 end
551 end
552 % -F: -sum{j~=i,k~=kstar} q(i,j,k,kstar,ni_c+2)*p2(i,ni_c+2,k,i,ni_c+2,k)
553 for j = 1:M
554 if j ~= i
555 for k2 = 1:K(i)
556 if k2 ~= kstar
557 qv = q_val(i, j, k2, kstar, ni_c+2);
558 idx = getP2Idx(i, ni_c+2, k2, i, ni_c+2, k2);
559 row(idx) = row(idx) - qv;
560 end
561 end
562 end
563 end
564 % -G: -sum{j~=i} q(i,j,kstar,kstar,ni_c+2)*p2(i,ni_c+2,kstar,i,ni_c+2,kstar)
565 for j = 1:M
566 if j ~= i
567 qv = q_val(i, j, kstar, kstar, ni_c+2);
568 idx = getP2Idx(i, ni_c+2, kstar, i, ni_c+2, kstar);
569 row(idx) = row(idx) - qv;
570 end
571 end
572 % -H: -sum{k~=kstar} q(i,i,k,kstar,ni_c+1)*p2(i,ni_c+1,k,i,ni_c+1,k)
573 for k2 = 1:K(i)
574 if k2 ~= kstar
575 qv = q_val(i, i, k2, kstar, ni_c+1);
576 idx = getP2Idx(i, ni_c+1, k2, i, ni_c+1, k2);
577 row(idx) = row(idx) - qv;
578 end
579 end
580 Aeq = [Aeq; row];
581 beq = [beq; 0];
582 end
583 end
584 end
585
586 %% COR1b: Correlation constraint boundary (from AMPL bnd_quadraticreduction_ld.mod)
587 % for each (i, kstar)
588 if verbose; fprintf(' COR1b constraints...\n'); end
589 for i = 1:M
590 for kstar = 1:K(i)
591 row = zeros(1, nVars);
592 % A': sum{j~=i,kj,hj,u~=kstar} q(j,i,kj,hj,1)*p2(j,1,kj,i,N-1,u)
593 for j = 1:M
594 if j ~= i
595 for kj = 1:K(j)
596 for hj = 1:K(j)
597 for u = 1:K(i)
598 if u ~= kstar
599 qv = q_val(j, i, kj, hj, 1);
600 idx = getP2Idx(j, 1, kj, i, N-1, u);
601 row(idx) = row(idx) + qv;
602 end
603 end
604 end
605 end
606 end
607 end
608 % C': sum{k~=kstar} q(i,i,kstar,k,N)*p2(i,N,kstar,i,N,kstar)
609 for k2 = 1:K(i)
610 if k2 ~= kstar
611 qv = q_val(i, i, kstar, k2, N);
612 idx = getP2Idx(i, N, kstar, i, N, kstar);
613 row(idx) = row(idx) + qv;
614 end
615 end
616 % -D': -sum{j~=i,k~=kstar} q(i,j,k,k,N)*p2(i,N,k,i,N,k)
617 for j = 1:M
618 if j ~= i
619 for k2 = 1:K(i)
620 if k2 ~= kstar
621 qv = q_val(i, j, k2, k2, N);
622 idx = getP2Idx(i, N, k2, i, N, k2);
623 row(idx) = row(idx) - qv;
624 end
625 end
626 end
627 end
628 % -E': -sum{j~=i,k~=kstar,h~=k} q(i,j,k,h,N)*p2(i,N,k,i,N,k)
629 for j = 1:M
630 if j ~= i
631 for k2 = 1:K(i)
632 if k2 ~= kstar
633 for h2 = 1:K(i)
634 if h2 ~= k2
635 qv = q_val(i, j, k2, h2, N);
636 idx = getP2Idx(i, N, k2, i, N, k2);
637 row(idx) = row(idx) - qv;
638 end
639 end
640 end
641 end
642 end
643 end
644 % -H': -sum{k~=kstar} q(i,i,k,kstar,N)*p2(i,N,k,i,N,k)
645 for k2 = 1:K(i)
646 if k2 ~= kstar
647 qv = q_val(i, i, k2, kstar, N);
648 idx = getP2Idx(i, N, k2, i, N, k2);
649 row(idx) = row(idx) - qv;
650 end
651 end
652 Aeq = [Aeq; row];
653 beq = [beq; 0];
654 end
655 end
656
657 %% THM4: QMIN inequality constraint
658 % for each (j, k, i):
659 % sum{t, h, nj, nt} nt*p2(j,nj,k,t,nt,h) - N*sum{h,nj,ni} p2(j,nj,k,i,ni,h) >= 0
660 if verbose; fprintf(' THM4 (QMIN) constraints...\n'); end
661 for j = 1:M
662 for kj = 1:K(j)
663 for i = 1:M
664 row = zeros(1, nVars);
665 % LHS: sum{t,h,nj,nt} nt*p2(j,nj,k,t,nt,h)
666 for t = 1:M
667 for ht = 1:K(t)
668 for nj = 0:N
669 for nt = 0:N
670 idx = getP2Idx(j, nj, kj, t, nt, ht);
671 row(idx) = row(idx) + nt;
672 end
673 end
674 end
675 end
676 % RHS: -N * sum{h,nj,ni} p2(j,nj,k,i,ni,h)
677 for hi = 1:K(i)
678 for nj = 0:N
679 for ni = 0:N
680 idx = getP2Idx(j, nj, kj, i, ni, hi);
681 row(idx) = row(idx) - N;
682 end
683 end
684 end
685 % >= 0 means -row*x <= 0
686 Aineq = [Aineq; -row];
687 bineq = [bineq; 0];
688 end
689 end
690 end
691
692 %% Build objective function
693 if verbose; fprintf('Building objective function...\n'); end
694 c = zeros(nVars, 1);
695
696 % Objective: p2(objective_queue, objective_n, objective_phase,
697 % objective_queue, objective_n, objective_phase)
698 idx_obj = getP2Idx(objective_queue, objective_n, objective_phase, ...
699 objective_queue, objective_n, objective_phase);
700 c(idx_obj) = 1;
701
702 if strcmp(sense, 'max')
703 c = -c; % linprog minimizes, so negate for maximization
704 end
705
706 %% Solve LP
707 if verbose
708 fprintf('Solving LP with %d variables and %d equality + %d inequality constraints...\n', ...
709 nVars, size(Aeq, 1), size(Aineq, 1));
710 end
711
712 if verbose
713 options = optimoptions('linprog', 'Display', 'final', 'Algorithm', 'interior-point');
714 else
715 options = optimoptions('linprog', 'Display', 'off', 'Algorithm', 'interior-point');
716 end
717
718 [x, fval, exitflag] = linprog(c, Aineq, bineq, Aeq, beq, lb, ub, options);
719
720 if strcmp(sense, 'max')
721 fval = -fval;
722 end
723
724 %% Extract results
725 result = struct();
726 result.objective = fval;
727 result.exitflag = exitflag;
728
729 if exitflag > 0
730 % Extract marginal probabilities p2(j,nj,k,j,nj,k)
731 result.p2marginals = zeros(M, N+1, max(K));
732 for j = 1:M
733 for nj = 0:N
734 for kj = 1:K(j)
735 idx = getP2Idx(j, nj, kj, j, nj, kj);
736 result.p2marginals(j, nj+1, kj) = x(idx);
737 end
738 end
739 end
740 end
741
742 if verbose; fprintf('\n=== Results ===\n'); end
743 if verbose; fprintf('Objective value: %f\n', fval); end
744 if verbose; fprintf('Exit flag: %d\n', exitflag); end
745 if exitflag > 0 && verbose
746 fprintf('\nMarginal probabilities p2(j,nj,k,j,nj,k):\n');
747 for j = 1:M
748 for nj = 0:N
749 for kj = 1:K(j)
750 val = result.p2marginals(j, nj+1, kj);
751 if val > 1e-8
752 fprintf(' p2(%d,%d,%d) = %.6f\n', j, nj, kj, val);
753 end
754 end
755 end
756 end
757 end
758end
759
760function qv = q_func(i, j, k, h, n, mu, v, alpha, r)
761% Q_FUNC - Compute load-dependent transition rate (1-based indices)
762%
763% i: source queue (1-based)
764% j: destination queue (1-based)
765% k: source phase (1-based)
766% h: destination phase (1-based)
767% n: population at queue i (0 returns 0)
768 if n == 0
769 qv = 0;
770 return;
771 end
772
773 % Load-dependent scaling factor
774 if n <= size(alpha, 2)
775 alpha_val = alpha(i, n);
776 else
777 alpha_val = 1;
778 end
779
780 if j ~= i
781 qv = r(i, j) * mu{i}(k, h) * alpha_val;
782 else
783 qv = (v{i}(k, h) + r(i, i) * mu{i}(k, h)) * alpha_val;
784 end
785end