LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
solver_mam_ag.m
1function [QN, UN, RN, TN, CN, XN, iter] = solver_mam_ag(sn, options)
2% SOLVER_MAM_AG AG methods for SolverMAM
3%
4% [QN, UN, RN, TN, CN, XN, ITER] = SOLVER_MAM_AG(SN, OPTIONS)
5%
6% Uses RCAT (Reversed Compound Agent Theorem) to find product-form
7% solutions for queueing networks.
8%
9% Methods:
10% 'inap' - Iterative Numerical Approximation Procedure (default, fast)
11% 'inapplus' - Improved INAP with weighted rates (no normalization)
12% 'exact' - Not available (autocat moved to line-legacy.git)
13%
14% Copyright (c) 2012-2025, Imperial College London
15% All rights reserved.
16
17M = sn.nstations;
18K = sn.nclasses;
19
20% Set default max states for truncation
21if isfield(options, 'config') && isfield(options.config, 'maxStates')
22 maxStates = options.config.maxStates;
23else
24 maxStates = 100;
25end
26
27% Set default tolerances
28if isfield(options, 'iter_tol') && ~isempty(options.iter_tol)
29 tol = options.iter_tol;
30else
31 tol = 1e-6;
32end
33
34if isfield(options, 'iter_max') && ~isempty(options.iter_max)
35 maxiter = options.iter_max;
36else
37 maxiter = 1000;
38end
39
40% Build RCAT model from network structure
41[R, AP, processMap, actionMap, N] = build_rcat(sn, maxStates);
42
43% Check if we have a valid model
44numProcesses = max(processMap(:));
45numActions = size(AP, 1);
46
47% Return early only if no processes found
48if numProcesses == 0
49 line_warning(mfilename, 'Network could not be mapped to RCAT format (no processes found).\n');
50 QN = zeros(M, K);
51 UN = zeros(M, K);
52 RN = zeros(M, K);
53 TN = zeros(M, K);
54 CN = zeros(1, K);
55 XN = zeros(1, K);
56 iter = 0;
57 return;
58end
59
60% If no actions but we have processes, solve using local rates only
61% This handles single-queue G-networks (Source -> Queue -> Sink)
62if numActions == 0
63 % No inter-station actions: solve equilibrium using only L matrices
64 x = [];
65 pi = cell(1, numProcesses);
66 Q = cell(1, numProcesses);
67 for p = 1:numProcesses
68 L = R{1, p}; % Local rate matrix is in R{numActions+1, p} = R{1, p} when numActions=0
69 % Convert to valid generator matrix
70 Qp = L - diag(L * ones(size(L, 1), 1));
71 Q{p} = ctmc_makeinfgen(Qp);
72 % Solve for equilibrium
73 pi{p} = ctmc_solve(Q{p});
74 end
75 iter = 0;
76 [QN, UN, RN, TN, CN, XN] = rcat_metrics(sn, x, pi, Q, processMap, actionMap, N);
77 return;
78end
79
80% Choose solver method
81method = options.method;
82if strcmp(method, 'default')
83 method = 'inap';
84end
85
86switch method
87 case 'inap'
88 % Fast iterative heuristic
89 [x, pi, Q, iter] = inap(R, AP, tol, maxiter, 'inap');
90
91 case 'inapplus'
92 % Improved INAP with weighted rates (no normalization)
93 [x, pi, Q, iter] = inap(R, AP, tol, maxiter, 'inapplus');
94
95 case 'exact'
96 % Optimization-based solver using autocat (not available in this version)
97 line_warning(mfilename, '''exact'' method not available. Falling back to inap.\n');
98 [x, pi, Q, iter] = inap(R, AP, tol, maxiter, 'inap');
99
100 otherwise
101 line_error(mfilename, 'Unknown method: %s\n', method);
102end
103
104% Convert RCAT solution to LINE metrics
105[QN, UN, RN, TN, CN, XN] = rcat_metrics(sn, x, pi, Q, processMap, actionMap, N);
106
107end
108
109%% Local Functions
110
111function [x, pi, Q, iter] = inap(R, AP, tol, maxiter, method)
112% INAP Iterative Numerical Approximation Procedure for RCAT
113%
114% Methods:
115% 'inap': x(a) = mean(Aa(i,j) * pi(i) / pi(j))
116% 'inapplus': x(a) = sum(Aa(i,j) * pi(i))
117
118if nargin < 5 || isempty(method)
119 method = 'inap';
120end
121
122% Parse R and AP
123A = size(AP, 1);
124ACT = AP(:, 1);
125PSV = AP(:, 2);
126numProcesses = max(AP(:));
127
128% Extract rate matrices
129Aa = cell(1, A);
130Pb = cell(1, A);
131for a = 1:A
132 Aa{a} = R{a, 1};
133 Pb{a} = R{a, 2};
134end
135
136% Extract local rates
137L = cell(1, numProcesses);
138for k = 1:numProcesses
139 L{k} = R{A+1, k};
140end
141
142% Get state space sizes
143N = zeros(1, numProcesses);
144for k = 1:numProcesses
145 N(k) = size(L{k}, 1);
146end
147
148% Initialize action rates randomly
149x = rand(A, 1);
150
151% Compute initial equilibrium
152[pi, Q] = compute_equilibrium(x, Aa, Pb, L, ACT, PSV, numProcesses, A, N);
153
154iter = 1;
155while iter <= maxiter
156 % Save old iteration result
157 piprev = pi;
158
159 % Update each action rate
160 for a = 1:A
161 k = ACT(a);
162
163 if strcmp(method, 'inapplus')
164 % inapplus: LAMBDA(i,j) = Aa{a}(i,j) * pi{k}(i)
165 % x(a) = sum(LAMBDA) for non-zero entries
166 LAMBDA_sum = 0;
167 for i = 1:N(k)
168 for j = 1:N(k)
169 if Aa{a}(i,j) > 0
170 LAMBDA_sum = LAMBDA_sum + Aa{a}(i,j) * pi{k}(i);
171 end
172 end
173 end
174 if LAMBDA_sum > 0
175 x(a) = LAMBDA_sum;
176 end
177 else
178 % inap: LAMBDA(i,j) = Aa{a}(i,j) * pi{k}(i) / pi{k}(j)
179 % x(a) = mean(LAMBDA) for non-zero entries
180 LAMBDA_vec = [];
181 for i = 1:N(k)
182 for j = 1:N(k)
183 if Aa{a}(i,j) > 0 && pi{k}(j) > 0
184 LAMBDA_vec(end+1) = Aa{a}(i,j) * pi{k}(i) / pi{k}(j);
185 end
186 end
187 end
188 if ~isempty(LAMBDA_vec)
189 x(a) = mean(LAMBDA_vec);
190 end
191 end
192 end
193
194 % Recompute equilibrium with new x
195 [pi, Q] = compute_equilibrium(x, Aa, Pb, L, ACT, PSV, numProcesses, A, N);
196
197 % Convergence check
198 maxErr = 0;
199 for k = 1:numProcesses
200 maxErr = max([maxErr, norm(pi{k} - piprev{k}, 1)]);
201 end
202
203 if maxErr < tol
204 return
205 end
206
207 iter = iter + 1;
208end
209
210end
211
212function [pi, Q] = compute_equilibrium(x, Aa, Pb, L, ACT, PSV, numProcesses, A, N)
213% Compute equilibrium distribution for each process given action rates x
214
215Q = cell(1, numProcesses);
216pi = cell(1, numProcesses);
217
218for k = 1:numProcesses
219 % Start with local/hidden rates
220 Qk = L{k} - diag(L{k} * ones(N(k), 1));
221
222 % Add contributions from each action
223 for c = 1:A
224 if PSV(c) == k
225 % Process k is passive for action c: add x(c) * Pb{c}
226 Qk = Qk + x(c) * Pb{c} - diag(Pb{c} * ones(N(k), 1));
227 elseif ACT(c) == k
228 % Process k is active for action c: add Aa{c}
229 Qk = Qk + Aa{c} - diag(Aa{c} * ones(N(k), 1));
230 end
231 end
232
233 % Convert to valid infinitesimal generator
234 Q{k} = ctmc_makeinfgen(Qk);
235
236 % Solve for equilibrium distribution using birth-death recursion for
237 % tridiagonal generators (numerically stable for large state spaces),
238 % falling back to ctmc_solve otherwise
239 if is_tridiagonal(Q{k})
240 pi{k} = birth_death_solve(Q{k});
241 else
242 pi{k} = ctmc_solve(Q{k});
243 end
244end
245
246end
247
248function [R, AP, processMap, actionMap, N] = build_rcat(sn, maxStates)
249% BUILD_RCAT Convert LINE network structure to RCAT format
250
251if nargin < 2
252 maxStates = 100;
253end
254
255M = sn.nstations;
256K = sn.nclasses;
257rt = sn.rt; % (M*K) x (M*K) routing table
258
259% Identify station types
260sourceStations = [];
261sinkStations = [];
262queueStations = [];
263
264for ist = 1:M
265 nodeIdx = sn.stationToNode(ist);
266 if sn.nodetype(nodeIdx) == NodeType.Source
267 sourceStations(end+1) = ist;
268 elseif sn.nodetype(nodeIdx) == NodeType.Sink
269 sinkStations(end+1) = ist;
270 else
271 % Queue, Delay, or other service stations
272 queueStations(end+1) = ist;
273 end
274end
275
276% Create process mapping: each (station, class) pair at queue stations
277% Note: Signal classes (negative customers) don't create separate processes
278% as they only modify the state of positive customer processes
279processIdx = 0;
280processMap = zeros(M, K);
281for ist = queueStations
282 for r = 1:K
283 % Skip Signal classes - they don't have their own queue state
284 if sn.issignal(r)
285 continue;
286 end
287 % Check if this station serves this class
288 if ~isnan(sn.rates(ist, r)) && sn.rates(ist, r) > 0
289 processIdx = processIdx + 1;
290 processMap(ist, r) = processIdx;
291 end
292 end
293end
294numProcesses = processIdx;
295
296if numProcesses == 0
297 R = {};
298 AP = [];
299 actionMap = [];
300 N = [];
301 return;
302end
303
304% Check for signal classes - G-networks with signals
305% G-network signals (negative customers, catastrophes, batch removal) are handled by:
306% 1. Adding signal arrival effects to the L matrix for positive customer processes
307% 2. Creating actions for signal routing between stations (job removal at destination)
308
309% Determine number of states for each process
310N = zeros(1, numProcesses);
311for p = 1:numProcesses
312 [ist, r] = find(processMap == p);
313 if ~isempty(ist)
314 ist = ist(1); r = r(1);
315 if sn.njobs(r) < Inf % Closed class
316 N(p) = sn.njobs(r) + 1; % States 0, 1, ..., njobs
317 else % Open class
318 N(p) = maxStates; % Truncate at maxStates
319 end
320 end
321end
322
323% Count actions: each routing transition (i,r) -> (j,s) where P > 0
324actionIdx = 0;
325actionMap = struct('from_station', {}, 'from_class', {}, ...
326 'to_station', {}, 'to_class', {}, 'prob', {}, ...
327 'isNegative', {}, 'isCatastrophe', {}, 'removalDistribution', {});
328
329for ist = queueStations
330 for r = 1:K
331 if processMap(ist, r) > 0
332 % Check if class r is a negative signal class
333 isNegativeClass = false;
334 isCatastropheClass = false;
335 removalDist = [];
336 if sn.issignal(r) && ~isnan(sn.signaltype{r}) && sn.signaltype{r} == SignalType.NEGATIVE
337 isNegativeClass = true;
338 % Check if class r is a catastrophe signal
339 if isfield(sn, 'isCatastrophe') && ~isempty(sn.isCatastrophe) && sn.isCatastrophe(r) > 0
340 isCatastropheClass = true;
341 end
342 % Get removal distribution for this class
343 if isfield(sn, 'signalRemovalDist') && ~isempty(sn.signalRemovalDist) && r <= length(sn.signalRemovalDist)
344 removalDist = sn.signalRemovalDist{r};
345 end
346 end
347
348 for jst = queueStations
349 for s = 1:K
350 if processMap(jst, s) > 0
351 % Get routing probability
352 prob_ij_rs = rt((ist-1)*K + r, (jst-1)*K + s);
353 if prob_ij_rs > 0 && (ist ~= jst || r ~= s)
354 % This is an action (departure from i,r triggers arrival at j,s)
355 actionIdx = actionIdx + 1;
356 actionMap(actionIdx).from_station = ist;
357 actionMap(actionIdx).from_class = r;
358 actionMap(actionIdx).to_station = jst;
359 actionMap(actionIdx).to_class = s;
360 actionMap(actionIdx).prob = prob_ij_rs;
361 actionMap(actionIdx).isNegative = isNegativeClass;
362 actionMap(actionIdx).isCatastrophe = isCatastropheClass;
363 actionMap(actionIdx).removalDistribution = removalDist;
364 end
365 end
366 end
367 end
368 end
369 end
370end
371numActions = actionIdx;
372
373% Initialize R and AP
374R = cell(numActions + 1, max(numProcesses, 2));
375if numActions > 0
376 AP = zeros(numActions, 2);
377else
378 AP = zeros(0, 2); % Empty matrix when no actions
379end
380
381% Identify sink nodes (nodetype = -1 = NodeType.Sink)
382% Use row vector to ensure for-loop doesn't execute when empty
383sinkNodes = find(sn.nodetype == NodeType.Sink)';
384if isempty(sinkNodes)
385 sinkNodes = []; % Ensure empty row vector, not column
386end
387
388% Build local/hidden rate matrices L for each process (R{end,k})
389for p = 1:numProcesses
390 [ist, r] = find(processMap == p);
391 ist = ist(1); r = r(1);
392 R{numActions + 1, p} = build_local_rates(sn, ist, r, N(p), rt, sourceStations, sinkNodes, K);
393end
394
395% Build active and passive matrices for each action
396for a = 1:numActions
397 am = actionMap(a);
398
399 % Active process (departure)
400 ist = am.from_station;
401 r = am.from_class;
402 p_active = processMap(ist, r);
403 AP(a, 1) = p_active;
404
405 % Get service rate
406 mu_ir = sn.rates(ist, r);
407 prob = am.prob;
408
409 % Active matrix: transition n -> n-1 with rate mu*prob (service completion)
410 Aa = zeros(N(p_active));
411 for n = 2:N(p_active)
412 Aa(n, n-1) = mu_ir * prob;
413 end
414 % Boundary: at max capacity
415 Aa(N(p_active), N(p_active)) = mu_ir * prob;
416 R{a, 1} = Aa;
417
418 % Passive process (arrival or signal effect)
419 jst = am.to_station;
420 s = am.to_class;
421 p_passive = processMap(jst, s);
422 AP(a, 2) = p_passive;
423
424 Pb = zeros(N(p_passive));
425 if am.isNegative
426 % NEGATIVE: Job removal at destination (G-network negative customer)
427 if am.isCatastrophe
428 % CATASTROPHE: All jobs are removed - all states transition to 1 (empty)
429 for n = 1:N(p_passive)
430 Pb(n, 1) = 1;
431 end
432 elseif ~isempty(am.removalDistribution)
433 % BATCH REMOVAL: Remove a random number of jobs based on distribution
434 % P[n, m] = probability of transition from n to m jobs
435 dist = am.removalDistribution;
436 for n = 1:N(p_passive)
437 if n == 1
438 % Empty queue: no effect
439 Pb(1, 1) = 1;
440 else
441 % For each possible resulting state m (from 1 to n)
442 for m = 1:n
443 k = n - m; % Number of jobs to remove to go from n to m
444 % Probability of removing exactly k jobs when queue has n-1 jobs (0-indexed)
445 if m > 1
446 % Remove exactly k jobs: P(removal = k)
447 prob = dist.evalPMF(k);
448 if prob > 0
449 Pb(n, m) = Pb(n, m) + prob;
450 end
451 else
452 % Remove all jobs (m = 1, i.e., state 0): P(removal >= n-1)
453 % = 1 - CDF(n-2) = 1 - sum_{j=0}^{n-2} P(removal = j)
454 cdfNMinus1 = 0;
455 for j = 0:(n-2)
456 cdfNMinus1 = cdfNMinus1 + dist.evalPMF(j);
457 end
458 probAtLeastN = 1 - cdfNMinus1;
459 if probAtLeastN > 0
460 Pb(n, 1) = Pb(n, 1) + probAtLeastN;
461 end
462 end
463 end
464 end
465 end
466 else
467 % DEFAULT: Remove exactly 1 job (original behavior)
468 % Empty queue: no effect (state 1 stays at state 1)
469 Pb(1, 1) = 1;
470 % Non-empty queues: decrement (n -> n-1)
471 for n = 2:(N(p_passive) - 1)
472 Pb(n, n-1) = 1;
473 end
474 % Boundary at max capacity: decrement
475 if N(p_passive) > 1
476 Pb(N(p_passive), N(p_passive) - 1) = 1;
477 end
478 end
479 else
480 % POSITIVE: Normal job arrival at destination
481 for n = 1:(N(p_passive) - 1)
482 Pb(n, n+1) = 1;
483 end
484 % Boundary: at max capacity
485 Pb(N(p_passive), N(p_passive)) = 1;
486 end
487 R{a, 2} = Pb;
488end
489
490end
491
492function L = build_local_rates(sn, ist, r, Np, rt, sourceStations, sinkNodes, K)
493% Build local/hidden transition matrix for process at station ist, class r
494% Note: sinkNodes contains node indices (not station indices) for Sink nodes
495
496L = zeros(Np);
497
498% External arrivals from source - separate positive, negative (single), batch, and catastrophe
499lambda_ir_pos = 0; % Positive arrivals
500lambda_ir_neg_single = 0; % Negative arrivals with single removal (default)
501lambda_ir_catastrophe = 0; % Catastrophe arrivals (remove all)
502% Batch removal arrivals: cell array of {rate, distribution} pairs
503batchArrivals = {};
504
505for isrc = sourceStations
506 for s_src = 1:K
507 % Check if source class s_src is a signal
508 isSignal = sn.issignal(s_src);
509
510 if isSignal
511 % For signals: they route to themselves (Signal -> Signal), but their effect
512 % is on positive customers at the destination station. We check if the signal
513 % routes to ANY class at this station (not just class r).
514 prob_src = 0;
515 for s_dst = 1:K
516 prob_src = prob_src + rt((isrc-1)*K + s_src, (ist-1)*K + s_dst);
517 end
518 else
519 % For regular classes: direct routing to (ist, r)
520 prob_src = rt((isrc-1)*K + s_src, (ist-1)*K + r);
521 end
522
523 if prob_src > 0 && ~isnan(sn.rates(isrc, s_src))
524 srcRate = sn.rates(isrc, s_src);
525 % Check if source class s_src is a negative or catastrophe signal
526 if isSignal && ~isnan(sn.signaltype{s_src}) && ...
527 (sn.signaltype{s_src} == SignalType.NEGATIVE || sn.signaltype{s_src} == SignalType.CATASTROPHE)
528 % Check if it's a catastrophe (either via isCatastrophe flag or signaltype)
529 isCat = (isfield(sn, 'isCatastrophe') && ~isempty(sn.isCatastrophe) && sn.isCatastrophe(s_src)) || ...
530 sn.signaltype{s_src} == SignalType.CATASTROPHE;
531 if isCat
532 lambda_ir_catastrophe = lambda_ir_catastrophe + srcRate * prob_src;
533 else
534 % Check if it has a removal distribution
535 removalDist = [];
536 if isfield(sn, 'signalRemovalDist') && ~isempty(sn.signalRemovalDist) && s_src <= length(sn.signalRemovalDist)
537 removalDist = sn.signalRemovalDist{s_src};
538 end
539 if ~isempty(removalDist)
540 batchArrivals{end+1} = {srcRate * prob_src, removalDist};
541 else
542 lambda_ir_neg_single = lambda_ir_neg_single + srcRate * prob_src;
543 end
544 end
545 else
546 lambda_ir_pos = lambda_ir_pos + srcRate * prob_src;
547 end
548 end
549 end
550end
551
552% Positive arrival transitions: n -> n+1
553if lambda_ir_pos > 0
554 for n = 1:(Np-1)
555 L(n, n+1) = lambda_ir_pos;
556 end
557end
558
559% Catastrophe arrival transitions: n -> 1 (for all n > 1)
560if lambda_ir_catastrophe > 0
561 for n = 2:Np
562 L(n, 1) = L(n, 1) + lambda_ir_catastrophe;
563 end
564end
565
566% Batch removal arrival transitions: n -> m at rate λ * P(remove n-m) for m < n
567for b = 1:length(batchArrivals)
568 rate = batchArrivals{b}{1};
569 dist = batchArrivals{b}{2};
570 for n = 2:Np
571 for m = 1:n
572 k = n - m; % Number of jobs to remove
573 if m > 1
574 % Remove exactly k jobs: P(removal = k)
575 prob = dist.evalPMF(k);
576 else
577 % Remove all jobs (m = 1): P(removal >= n-1)
578 cdfNMinus1 = 0;
579 for j = 0:(n-2)
580 cdfNMinus1 = cdfNMinus1 + dist.evalPMF(j);
581 end
582 prob = 1 - cdfNMinus1;
583 end
584 if prob > 0
585 L(n, m) = L(n, m) + rate * prob;
586 end
587 end
588 end
589end
590
591% Single removal negative arrival transitions: n -> n-1 (only if queue non-empty)
592if lambda_ir_neg_single > 0
593 for n = 2:Np
594 L(n, n-1) = L(n, n-1) + lambda_ir_neg_single;
595 end
596end
597
598% Service rate at this station
599mu_ir = sn.rates(ist, r);
600if ~isnan(mu_ir) && mu_ir > 0
601 % Departures to sink (use rtnodes with node indices)
602 % Get the node index for this station
603 nodeIdx = sn.stationToNode(ist);
604
605 prob_sink = 0;
606 if isfield(sn, 'rtnodes') && ~isempty(sn.rtnodes)
607 nNodes = length(sn.nodetype);
608 for jsnk = sinkNodes
609 for s = 1:K
610 % rtnodes indices: (nodeIdx-1)*K + classIdx
611 fromIdx = (nodeIdx - 1) * K + r;
612 toIdx = (jsnk - 1) * K + s;
613 if fromIdx <= size(sn.rtnodes, 1) && toIdx <= size(sn.rtnodes, 2)
614 prob_sink = prob_sink + sn.rtnodes(fromIdx, toIdx);
615 end
616 end
617 end
618 end
619
620 % Self-routing (stays at same station, same class)
621 prob_self = rt((ist-1)*K + r, (ist-1)*K + r);
622
623 % Combined local departure rate
624 local_departure_rate = mu_ir * prob_sink;
625
626 % Departure transitions: n -> n-1 (add to existing negative arrival effects)
627 for n = 2:Np
628 L(n, n-1) = L(n, n-1) + local_departure_rate;
629 end
630
631 % Self-service transitions: n -> n (diagonal, for phase transitions)
632 for n = 2:Np
633 L(n, n) = mu_ir * prob_self;
634 end
635end
636
637end
638
639function [QN, UN, RN, TN, CN, XN] = rcat_metrics(sn, x, pi, Q, processMap, actionMap, N)
640% RCAT_METRICS Convert RCAT solution to LINE performance metrics
641
642M = sn.nstations;
643K = sn.nclasses;
644
645QN = zeros(M, K);
646UN = zeros(M, K);
647RN = zeros(M, K);
648TN = zeros(M, K);
649
650% Compute metrics for each (station, class) pair
651for ist = 1:M
652 for r = 1:K
653 p = processMap(ist, r);
654 if p > 0 && ~isempty(pi) && p <= length(pi) && ~isempty(pi{p})
655 Np = length(pi{p});
656
657 % Queue length: E[N] = sum_{n=0}^{Np-1} n * pi(n+1)
658 QN(ist, r) = (0:(Np-1)) * pi{p}(:);
659
660 % Utilization: P(N > 0) = 1 - pi(0) = 1 - pi{p}(1)
661 UN(ist, r) = 1 - pi{p}(1);
662
663 % Throughput: compute from utilization and service rate
664 mu_ir = sn.rates(ist, r);
665 if ~isnan(mu_ir) && mu_ir > 0
666 TN(ist, r) = mu_ir * UN(ist, r);
667 end
668 end
669 end
670end
671
672% Handle self-looping classes: they always stay at their reference station
673% and share the server with other classes under PS scheduling
674if isfield(sn, 'isslc') && any(sn.isslc)
675 for r = 1:K
676 if sn.isslc(r)
677 refst = sn.refstat(r);
678 if refst > 0 && refst <= M
679 % Self-looping class: all jobs stay at reference station
680 QN(refst, r) = sn.njobs(r);
681
682 % Service rate for this class
683 mu_ir = sn.rates(refst, r);
684 if ~isnan(mu_ir) && mu_ir > 0
685 nservers = sn.nservers(refst);
686 if isinf(nservers)
687 nservers = 1; % Delay station
688 end
689
690 % For PS scheduling, compute utilization using offered load
691 % First get utilization from other classes at this station
692 other_util = 0;
693 for s = 1:K
694 if s ~= r && ~sn.isslc(s)
695 other_util = other_util + UN(refst, s);
696 end
697 end
698
699 % Remaining capacity is shared with SLC
700 % Total utilization should sum to offered load (capped at 1)
701 % For SLC with N jobs always present: U_slc = T_slc / mu
702 % Under PS, T_slc = mu * N_slc * (1 - fraction_to_other) / N_slc
703 % Simplified: remaining capacity goes to SLC
704 remaining_capacity = max(0, 1 - other_util);
705
706 % SLC utilization: min(demand, remaining capacity)
707 % Demand for SLC = N_slc / mu = QN / mu
708 slc_demand = QN(refst, r) / mu_ir;
709 UN(refst, r) = min(slc_demand, remaining_capacity);
710 TN(refst, r) = mu_ir * UN(refst, r);
711 end
712 end
713 end
714 end
715end
716
717% Response times from Little's law: R = Q / T
718for ist = 1:M
719 for r = 1:K
720 if TN(ist, r) > 0
721 RN(ist, r) = QN(ist, r) / TN(ist, r);
722 else
723 RN(ist, r) = 0;
724 end
725 end
726end
727
728% System metrics
729CN = zeros(1, K);
730XN = zeros(1, K);
731
732% For open classes: system throughput = arrival rate, system response time = sum of response times
733for r = 1:K
734 if sn.njobs(r) >= Inf % Open class
735 % System throughput equals arrival rate (from source)
736 for ist = 1:M
737 nodeIdx = sn.stationToNode(ist);
738 if sn.nodetype(nodeIdx) == NodeType.Source
739 XN(r) = sn.rates(ist, r);
740 break;
741 end
742 end
743 % System response time = sum over all stations
744 CN(r) = sum(RN(:, r));
745 else
746 % Closed class: use reference station
747 refst = sn.refstat(r);
748 if refst > 0 && refst <= M
749 XN(r) = TN(refst, r);
750 if XN(r) > 0
751 CN(r) = sn.njobs(r) / XN(r);
752 end
753 end
754 end
755end
756
757end
758
759function pi = birth_death_solve(Q)
760% BIRTH_DEATH_SOLVE Solve equilibrium of a birth-death (tridiagonal) CTMC
761%
762% For a birth-death chain with birth rate lambda_n = Q(n, n+1) and
763% death rate mu_n = Q(n, n-1), the equilibrium is computed using the
764% recursion pi(n) = pi(n-1) * lambda(n-1) / mu(n).
765%
766% This is numerically stable and avoids the ill-conditioned linear system
767% that plagues null-space methods for large state spaces.
768
769n = size(Q, 1);
770if n <= 1
771 pi = 1;
772 return;
773end
774
775pi = zeros(1, n);
776pi(1) = 1.0;
777
778for i = 2:n
779 birth_rate = Q(i-1, i);
780 death_rate = Q(i, i-1);
781 if death_rate > 0
782 pi(i) = pi(i-1) * birth_rate / death_rate;
783 else
784 pi(i) = 0;
785 end
786end
787
788total = sum(pi);
789if total > 0
790 pi = pi / total;
791else
792 pi = ones(1, n) / n;
793end
794
795end
796
797function result = is_tridiagonal(Q)
798% IS_TRIDIAGONAL Check if a matrix is tridiagonal
799n = size(Q, 1);
800result = true;
801for i = 1:n
802 for j = 1:n
803 if abs(i - j) > 1 && abs(Q(i, j)) > 1e-14
804 result = false;
805 return;
806 end
807 end
808end
809end
Definition mmt.m:92