1function [QN, UN, RN, TN, CN, XN, iter] = solver_mam_ag(sn, options)
2% SOLVER_MAM_AG AG methods
for SolverMAM
4% [QN, UN, RN, TN, CN, XN, ITER] = SOLVER_MAM_AG(SN, OPTIONS)
6% Uses RCAT (Reversed Compound Agent Theorem) to find product-form
7% solutions
for queueing networks.
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)
14% Copyright (c) 2012-2025, Imperial College London
20% Set
default max states
for truncation
21if isfield(options,
'config') && isfield(options.config,
'maxStates')
22 maxStates = options.config.maxStates;
27% Set default tolerances
28if isfield(options, 'iter_tol') && ~isempty(options.iter_tol)
29 tol = options.iter_tol;
34if isfield(options, 'iter_max') && ~isempty(options.iter_max)
35 maxiter = options.iter_max;
40% Build RCAT model from network structure
41[R, AP, processMap, actionMap, N] = build_rcat(sn, maxStates);
43% Check if we have a valid model
44numProcesses = max(processMap(:));
45numActions = size(AP, 1);
47% Return early only if no processes found
49 line_warning(mfilename, 'Network could not be mapped to RCAT format (no processes found).\n');
60% If no actions but we have processes, solve using local rates only
61% This handles single-queue G-networks (Source -> Queue -> Sink)
63 % No inter-station actions: solve equilibrium using only L matrices
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});
76 [QN, UN, RN, TN, CN, XN] = rcat_metrics(sn, x, pi, Q, processMap, actionMap, N);
81method = options.method;
82if strcmp(method,
'default')
88 % Fast iterative heuristic
89 [x, pi, Q, iter] = inap(R, AP, tol, maxiter, 'inap');
92 % Improved INAP with weighted rates (no normalization)
93 [x, pi, Q, iter] = inap(R, AP, tol, maxiter, 'inapplus');
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');
101 line_error(mfilename, 'Unknown method: %s\n', method);
104% Convert RCAT solution to LINE metrics
105[QN, UN, RN, TN, CN, XN] = rcat_metrics(sn, x, pi, Q, processMap, actionMap, N);
111function [x, pi, Q, iter] = inap(R, AP, tol, maxiter, method)
112% INAP Iterative Numerical Approximation Procedure for RCAT
115% 'inap': x(a) = mean(Aa(i,j) * pi(i) / pi(j))
116% 'inapplus': x(a) = sum(Aa(i,j) * pi(i))
118if nargin < 5 || isempty(method)
126numProcesses = max(AP(:));
128% Extract rate matrices
137L = cell(1, numProcesses);
138for k = 1:numProcesses
142% Get state space sizes
143N = zeros(1, numProcesses);
144for k = 1:numProcesses
145 N(k) = size(L{k}, 1);
148% Initialize action rates randomly
151% Compute initial equilibrium
152[pi, Q] = compute_equilibrium(x, Aa, Pb, L, ACT, PSV, numProcesses, A, N);
156 % Save old iteration result
159 % Update each action rate
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
170 LAMBDA_sum = LAMBDA_sum + Aa{a}(i,j) * pi{k}(i);
178 % inap: LAMBDA(i,j) = Aa{a}(i,j) * pi{k}(i) / pi{k}(j)
179 % x(a) = mean(LAMBDA)
for non-zero entries
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);
188 if ~isempty(LAMBDA_vec)
189 x(a) = mean(LAMBDA_vec);
194 % Recompute equilibrium with
new x
195 [pi, Q] = compute_equilibrium(x, Aa, Pb, L, ACT, PSV, numProcesses, A, N);
199 for k = 1:numProcesses
200 maxErr = max([maxErr, norm(pi{k} - piprev{k}, 1)]);
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
215Q = cell(1, numProcesses);
216pi = cell(1, numProcesses);
218for k = 1:numProcesses
219 % Start with local/hidden rates
220 Qk = L{k} - diag(L{k} * ones(N(k), 1));
222 % Add contributions from each action
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));
228 % Process k
is active
for action c: add Aa{c}
229 Qk = Qk + Aa{c} - diag(Aa{c} * ones(N(k), 1));
233 % Convert to valid infinitesimal generator
234 Q{k} = ctmc_makeinfgen(Qk);
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});
242 pi{k} = ctmc_solve(Q{k});
248function [R, AP, processMap, actionMap, N] = build_rcat(sn, maxStates)
249% BUILD_RCAT Convert LINE network structure to RCAT format
257rt = sn.rt; % (M*K) x (M*K) routing table
259% Identify station types
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;
271 % Queue, Delay, or other service stations
272 queueStations(end+1) = ist;
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
280processMap = zeros(M, K);
281for ist = queueStations
283 % Skip Signal
classes - they don't have their own queue state
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;
294numProcesses = processIdx;
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)
309% Determine number of states
for each process
310N = zeros(1, numProcesses);
311for p = 1:numProcesses
312 [ist, r] = find(processMap == p);
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
318 N(p) = maxStates; % Truncate at maxStates
323% Count actions: each routing transition (i,r) -> (j,s) where
P > 0
325actionMap =
struct(
'from_station', {},
'from_class', {}, ...
326 'to_station', {},
'to_class', {},
'prob', {}, ...
327 'isNegative', {},
'isCatastrophe', {},
'removalDistribution', {});
329for ist = queueStations
331 if processMap(ist, r) > 0
332 % Check
if class r
is a negative signal class
333 isNegativeClass =
false;
334 isCatastropheClass =
false;
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;
342 % Get removal distribution
for this class
343 if isfield(sn,
'signalRemovalDist') && ~isempty(sn.signalRemovalDist) && r <= length(sn.signalRemovalDist)
344 removalDist = sn.signalRemovalDist{r};
348 for jst = queueStations
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;
371numActions = actionIdx;
374R = cell(numActions + 1, max(numProcesses, 2));
376 AP = zeros(numActions, 2);
378 AP = zeros(0, 2); % Empty matrix when no actions
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)';
385 sinkNodes = []; % Ensure empty row vector, not column
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);
395% Build active and passive matrices
for each action
399 % Active process (departure)
400 ist = am.from_station;
402 p_active = processMap(ist, r);
406 mu_ir = sn.rates(ist, r);
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;
414 % Boundary: at max capacity
415 Aa(N(p_active), N(p_active)) = mu_ir * prob;
418 % Passive process (arrival or signal effect)
421 p_passive = processMap(jst, s);
422 AP(a, 2) = p_passive;
424 Pb = zeros(N(p_passive));
426 % NEGATIVE: Job removal at destination (G-network negative customer)
428 % CATASTROPHE: All jobs are removed - all states transition to 1 (empty)
429 for n = 1:N(p_passive)
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)
438 % Empty queue: no effect
441 % For each possible resulting state m (from 1 to 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)
446 % Remove exactly k jobs:
P(removal = k)
447 prob = dist.evalPMF(k);
449 Pb(n, m) = Pb(n, m) + prob;
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)
456 cdfNMinus1 = cdfNMinus1 + dist.evalPMF(j);
458 probAtLeastN = 1 - cdfNMinus1;
460 Pb(n, 1) = Pb(n, 1) + probAtLeastN;
467 % DEFAULT: Remove exactly 1 job (original behavior)
468 % Empty queue: no effect (state 1 stays at state 1)
470 % Non-empty queues: decrement (n -> n-1)
471 for n = 2:(N(p_passive) - 1)
474 % Boundary at max capacity: decrement
476 Pb(N(p_passive), N(p_passive) - 1) = 1;
480 % POSITIVE: Normal job arrival at destination
481 for n = 1:(N(p_passive) - 1)
484 % Boundary: at max capacity
485 Pb(N(p_passive), N(p_passive)) = 1;
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
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
505for isrc = sourceStations
507 % Check
if source
class s_src
is a signal
508 isSignal = sn.issignal(s_src);
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).
516 prob_src = prob_src + rt((isrc-1)*K + s_src, (ist-1)*K + s_dst);
519 % For regular
classes: direct routing to (ist, r)
520 prob_src = rt((isrc-1)*K + s_src, (ist-1)*K + r);
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;
532 lambda_ir_catastrophe = lambda_ir_catastrophe + srcRate * prob_src;
534 % Check if it has a removal distribution
536 if isfield(sn, 'signalRemovalDist
') && ~isempty(sn.signalRemovalDist) && s_src <= length(sn.signalRemovalDist)
537 removalDist = sn.signalRemovalDist{s_src};
539 if ~isempty(removalDist)
540 batchArrivals{end+1} = {srcRate * prob_src, removalDist};
542 lambda_ir_neg_single = lambda_ir_neg_single + srcRate * prob_src;
546 lambda_ir_pos = lambda_ir_pos + srcRate * prob_src;
552% Positive arrival transitions: n -> n+1
555 L(n, n+1) = lambda_ir_pos;
559% Catastrophe arrival transitions: n -> 1 (for all n > 1)
560if lambda_ir_catastrophe > 0
562 L(n, 1) = L(n, 1) + lambda_ir_catastrophe;
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};
572 k = n - m; % Number of jobs to remove
574 % Remove exactly k jobs: P(removal = k)
575 prob = dist.evalPMF(k);
577 % Remove all jobs (m = 1): P(removal >= n-1)
580 cdfNMinus1 = cdfNMinus1 + dist.evalPMF(j);
582 prob = 1 - cdfNMinus1;
585 L(n, m) = L(n, m) + rate * prob;
591% Single removal negative arrival transitions: n -> n-1 (only if queue non-empty)
592if lambda_ir_neg_single > 0
594 L(n, n-1) = L(n, n-1) + lambda_ir_neg_single;
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);
606 if isfield(sn, 'rtnodes
') && ~isempty(sn.rtnodes)
607 nNodes = length(sn.nodetype);
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);
620 % Self-routing (stays at same station, same class)
621 prob_self = rt((ist-1)*K + r, (ist-1)*K + r);
623 % Combined local departure rate
624 local_departure_rate = mu_ir * prob_sink;
626 % Departure transitions: n -> n-1 (add to existing negative arrival effects)
628 L(n, n-1) = L(n, n-1) + local_departure_rate;
631 % Self-service transitions: n -> n (diagonal, for phase transitions)
633 L(n, n) = mu_ir * prob_self;
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
650% Compute metrics for each (station, class) pair
653 p = processMap(ist, r);
654 if p > 0 && ~isempty(pi) && p <= length(pi) && ~isempty(pi{p})
657 % Queue length: E[N] = sum_{n=0}^{Np-1} n * pi(n+1)
658 QN(ist, r) = (0:(Np-1)) * pi{p}(:);
660 % Utilization: P(N > 0) = 1 - pi(0) = 1 - pi{p}(1)
661 UN(ist, r) = 1 - pi{p}(1);
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);
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)
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);
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);
687 nservers = 1; % Delay station
690 % For PS scheduling, compute utilization using offered load
691 % First get utilization from other classes at this station
694 if s ~= r && ~sn.isslc(s)
695 other_util = other_util + UN(refst, s);
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);
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);
717% Response times from Little's law: R = Q / T
721 RN(ist, r) = QN(ist, r) / TN(ist, r);
732% For open
classes: system throughput = arrival rate, system response time = sum of response times
734 if sn.njobs(r) >= Inf % Open
class
735 % System throughput equals arrival rate (from source)
737 nodeIdx = sn.stationToNode(ist);
738 if sn.nodetype(nodeIdx) == NodeType.Source
739 XN(r) = sn.rates(ist, r);
743 % System response time = sum over all stations
744 CN(r) = sum(RN(:, r));
746 % Closed
class: use reference station
747 refst = sn.refstat(r);
748 if refst > 0 && refst <= M
749 XN(r) = TN(refst, r);
751 CN(r) = sn.njobs(r) / XN(r);
759function pi = birth_death_solve(Q)
760% BIRTH_DEATH_SOLVE Solve equilibrium of a birth-death (tridiagonal) CTMC
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).
766% This
is numerically stable and avoids the ill-conditioned linear system
767% that plagues null-space methods for large state spaces.
779 birth_rate = Q(i-1, i);
780 death_rate = Q(i, i-1);
782 pi(i) = pi(i-1) * birth_rate / death_rate;
797function result = is_tridiagonal(Q)
798% IS_TRIDIAGONAL Check if a matrix
is tridiagonal
803 if abs(i - j) > 1 && abs(Q(i, j)) > 1e-14