LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
solver_mam_basic.m
1function [QN,UN,RN,TN,CN,XN,totiter] = solver_mam_basic(sn, options)
2% [Q,U,R,T,C,X] = SOLVER_MAM_BASIC(QN, PH, OPTIONS)
3
4% This solver uses MAM to solve queues in isolation, but simplifies the
5% traffic equations by using visits to rescale the flows into the queue
6% inputs.
7%
8% Copyright (c) 2012-2026, Imperial College London
9% All rights reserved.
10
11config = options.config;
12tol = options.tol;
13
14PH = sn.proc;
15%% generate local state spaces
16I = sn.nnodes;
17M = sn.nstations;
18K = sn.nclasses;
19C = sn.nchains;
20N = sn.njobs';
21V = cellsum(sn.visits);
22S = 1./sn.rates;
23Lchain = sn_get_demands_chain(sn);
24
25QN = zeros(M,K);
26UN = zeros(M,K);
27RN = zeros(M,K);
28TN = zeros(M,K);
29CN = zeros(1,K);
30XN = zeros(1,K);
31
32% Track stations using exact MAP/D/c solver (skip post-processing for these)
33mapdcStations = false(M, 1);
34
35pie = {};
36D0 = {};
37
38lambda = zeros(1,C);
39chainSysArrivals = cell(1,C);
40TN_1 = TN+Inf;
41
42it = 0;
43
44% open queueing system (one node is the external world)
45% first build the joint arrival process
46for ist=1:M
47 switch sn.sched(ist)
48 case {SchedStrategy.FCFS, SchedStrategy.HOL, SchedStrategy.FCFSPRPRIO, SchedStrategy.PS}
49 for k=1:K
50 % Skip Det processes - they will be handled by qsys_mapdc
51 if sn.procid(ist, k) == ProcessType.DET
52 % For Det, we don't need MAP representation since we use exact solver
53 pie{ist}{k} = [];
54 D0{ist,k} = [];
55 continue;
56 end
57 % divide service time by number of servers and put
58 % later a surrogate delay server in tandem to compensate
59 PH{ist}{k} = map_scale(PH{ist}{k}, S(ist,k)/sn.nservers(ist));
60 pie{ist}{k} = map_pie(PH{ist}{k});
61 D0{ist,k} = PH{ist}{k}{1};
62 if any(isnan(D0{ist,k}))
63 D0{ist,k} = -GlobalConstants.Immediate;
64 pie{ist}{k} = 1;
65 PH{ist}{k} = map_exponential(GlobalConstants.Immediate);
66 end
67 end
68 end
69end % i
70
71isOpen = false;
72isClosed = false;
73if any(isinf(sn.njobs))
74 isOpen = true;
75end
76if any(isfinite(sn.njobs))
77 isClosed = true;
78end
79isMixed = isOpen & isClosed;
80
81lambdas_inchain = cell(1,C);
82for c=1:C
83 inchain = sn.inchain{c};
84 lambdas_inchain{c} = sn.rates(sn.refstat(inchain(1)),inchain);
85 %lambdas_inchain{c} = lambdas_inchain{c}(isfinite(lambdas_inchain{c}));
86 lambda(c) = sum(lambdas_inchain{c}(isfinite(lambdas_inchain{c})));
87 ist = sn.refstat(inchain(1)); % identical for all classes in the chain
88 if isinf(sum(sn.njobs(inchain))) % if open chain
89 % ist here is the source
90 % assemble a MMAP for the arrival process from all classes
91 for k=1:K
92 if isnan(PH{ist}{k}{1})
93 PH{ist}{k} = map_exponential(Inf); % no arrivals from this class
94 end
95 end
96 inchain = sn.inchain{c};
97 k = inchain(1);
98 chainSysArrivals{c} = {PH{ist}{k}{1},PH{ist}{k}{2},PH{ist}{k}{2}};
99 for ki=2:length(inchain)
100 k = inchain(ki);
101 if isnan(PH{ist}{k}{1})
102 PH{ist}{k} = map_exponential(Inf); % no arrivals from this class
103 end
104 chainSysArrivals{c} = mmap_super_safe({chainSysArrivals{c},{PH{ist}{k}{1},PH{ist}{k}{2},PH{ist}{k}{2}}}, config.space_max, 'default');
105 end
106 TN(ist,inchain') = lambdas_inchain{c};
107 end
108end
109
110sd = isfinite(sn.nservers);
111
112while max(max(abs(TN-TN_1))) > tol && it <= options.iter_max %#ok<max>
113 it = it + 1;
114 TN_1 = TN;
115 Umax = max(sum(UN(sd,:),2));
116 if Umax >=1
117 lambda = lambda * 1/Umax;
118 else
119 for c=1:C
120 inchain = sn.inchain{c};
121 if ~isinf(sum(sn.njobs(inchain))) % if closed chain
122 Nc = sum(sn.njobs(inchain)); % closed population
123 QNc = max(tol, sum(sum(QN(:,inchain),2,"omitnan"))); %#ok<NANSUM>
124 TNlb = Nc./sum(Lchain(:,c));
125 if it == 1
126 lambda(c) = TNlb; % lower bound
127 else
128 lambda(c) = lambda(c) * it/options.iter_max + (Nc / QNc) * lambda(c) * (options.iter_max-it)/options.iter_max; % iteration-averaged regula falsi;
129 end
130 end
131 end
132 end
133
134 for c=1:C
135 inchain = sn.inchain{c};
136 chainSysArrivals{c} = mmap_exponential(lambda(c));
137 for m=1:M
138 TN(m,inchain) = V(m,inchain) .* lambda(c);
139 end
140 end
141
142 for ind=1:I
143 if sn.isstation(ind)
144 ist = sn.nodeToStation(ind);
145 switch sn.nodetype(ind)
146 case NodeType.Join
147 for c=1:C
148 inchain = sn.inchain{c};
149 for k=inchain
150 fanin = nnz(sn.rtnodes(:, (ind-1)*K+k));
151 TN(ist,k) = lambda(c)*V(ist,k)/fanin;
152 UN(ist,k) = 0;
153 QN(ist,k) = 0;
154 RN(ist,k) = 0;
155 end
156 end
157 otherwise
158 switch sn.sched(ist)
159 case SchedStrategy.INF
160 for c=1:C
161 inchain = sn.inchain{c};
162 for k=inchain
163 TN(ist,k) = lambda(c)*V(ist,k);
164 UN(ist,k) = S(ist,k)*TN(ist,k);
165 QN(ist,k) = TN(ist,k).*S(ist,k)*V(ist,k);
166 RN(ist,k) = QN(ist,k)/TN(ist,k);
167 end
168 end
169 case SchedStrategy.PS
170 for c=1:C
171 inchain = sn.inchain{c};
172 for k=inchain
173 TN(ist,k) = lambda(c)*V(ist,k);
174 UN(ist,k) = S(ist,k)*TN(ist,k);
175 end
176 %Nc = sum(sn.njobs(inchain)); % closed population
177 Uden = min([1-GlobalConstants.FineTol,sum(UN(ist,:))]);
178 for k=inchain
179 %QN(ist,k) = (UN(ist,k)-UN(ist,k)^(Nc+1))/(1-Uden); % geometric bound type approximation
180 QN(ist,k) = UN(ist,k)/(1-Uden);
181 RN(ist,k) = QN(ist,k)/TN(ist,k);
182 end
183 end
184 case {SchedStrategy.FCFS, SchedStrategy.HOL, SchedStrategy.FCFSPRPRIO}
185 chainArrivalAtNode = cell(1,C);
186 rates = cell(M,C);
187 mapdcUsed = false; % Flag for exact MAP/D/c solver
188 for c=1:C %for each chain
189 rates{ist,c} = V(ist,:) .* lambda(c); % visits of classes within the chain
190 inchain = find(sn.chains(c,:))';
191 markProb = rates{ist,c}(inchain) / sum(rates{ist,c}(inchain));
192 markProb(isnan(markProb)) = 0;
193 chainArrivalAtNode{c} = mmap_mark(chainSysArrivals{c}, markProb);
194 chainArrivalAtNode{c} = mmap_normalize(chainArrivalAtNode{c});
195 chainArrivalAtNode{c} = mmap_scale(chainArrivalAtNode{c}, 1./rates{ist,c}, 0); % non-iterative approximation
196 if c == 1
197 aggrArrivalAtNode = mmap_super_safe({chainArrivalAtNode{c}, mmap_exponential(0,1)}, config.space_max, 'default');
198 aggrArrivalAtNode = {aggrArrivalAtNode{1} aggrArrivalAtNode{2} aggrArrivalAtNode{2}};
199 lc = map_lambda(chainArrivalAtNode{c});
200 if lc>0
201 aggrArrivalAtNode = mmap_scale(aggrArrivalAtNode, 1/lc, 0); % non-iterative approximation
202 end
203 else
204 aggrArrivalAtNode = mmap_super_safe({aggrArrivalAtNode, chainArrivalAtNode{c}}, config.space_max, 'default');
205 end
206 end
207 Qret = cell(1,K);
208 if (strcmp(sn.sched(ist),SchedStrategy.HOL) && any(sn.classprio ~= sn.classprio(1))) % if priorities are not identical
209 [uK,iK] = unique(sn.classprio);
210 % BUTools convention: D1=lowest priority, DK=highest priority
211 % LINE convention: lower value = higher priority
212 % unique() returns ascending order, so we need to reverse for BUTools
213 iK = flipud(iK(:));
214 if length(uK) == length(sn.classprio) % if all priorities are different
215 [Qret{iK'}] = MMAPPH1NPPR({aggrArrivalAtNode{[1;2+iK]}}, {pie{ist}{iK}}, {D0{ist,iK}}, 'ncMoms', 1);
216 else
217 line_error(mfilename,'Solver MAM requires either identical priorities or all distinct priorities');
218 end
219 elseif (sn.sched(ist)==SchedStrategy.FCFSPRPRIO && any(sn.classprio ~= sn.classprio(1))) % FCFS preemptive resume priority
220 [uK,iK] = unique(sn.classprio);
221 % BUTools convention: D1=lowest priority, DK=highest priority
222 % LINE convention: lower value = higher priority
223 % unique() returns ascending order, so we need to reverse for BUTools
224 iK = flipud(iK(:));
225 if length(uK) == length(sn.classprio) % if all priorities are different
226 [Qret{iK'}] = MMAPPH1PRPR({aggrArrivalAtNode{[1;2+iK]}}, {pie{ist}{iK}}, {D0{ist,iK}}, 'ncMoms', 1);
227 else
228 line_error(mfilename,'Solver MAM requires either identical priorities or all distinct priorities');
229 end
230 else
231 aggrUtil = sum(mmap_lambda(aggrArrivalAtNode)./(GlobalConstants.FineTol+sn.rates(ist,1:K)*sn.nservers(ist)), 'omitnan'); %% to debug
232 aggrLambda = mmap_lambda(aggrArrivalAtNode);
233 if aggrUtil < 1-GlobalConstants.FineTol
234 if any(isinf(N))
235 % Check for MAP/D/c: single-class open model with deterministic service
236 isMapDc = (K == 1) && (sn.procid(ist, 1) == ProcessType.DET);
237 if isMapDc
238 % Use exact MAP/D/c solver from Q-MAM
239 D0_arr = aggrArrivalAtNode{1};
240 D1_arr = aggrArrivalAtNode{2};
241 detServiceTime = S(ist, 1); % 1/rate = service time
242 numServers = sn.nservers(ist);
243
244 mapdcResult = qsys_mapdc(D0_arr, D1_arr, detServiceTime, numServers);
245 Qret{1} = mapdcResult.meanQueueLength;
246 % Store result for later use to skip surrogate delay adjustment
247 mapdcUsed = true;
248 mapdcStations(ist) = true; % Mark for skipping post-processing
249 line_debug('Using exact MAP/D/%d solver: Q=%.4f, W=%.4f', ...
250 numServers, mapdcResult.meanQueueLength, mapdcResult.meanWaitingTime);
251 else
252 mapdcUsed = false;
253 [Qret{1:K}] = MMAPPH1FCFS({aggrArrivalAtNode{[1,3:end]}}, {pie{ist}{:}}, {D0{ist,:}}, 'ncMoms', 1);
254 end
255 else % all closed classes
256 maxLevel = sum(N(isfinite(N)))+1;
257 D = {aggrArrivalAtNode{[1,3:end]}};
258 pdistr_k = cell(1,K);
259 if map_lambda(D)< GlobalConstants.FineTol
260 for k=1:K
261 pdistr_k = [1-GlobalConstants.FineTol, GlobalConstants.FineTol];
262 Qret{k} = GlobalConstants.FineTol / sn.rates(ist);
263 end
264 else
265 if sn.isfunction(ist)
266 alpharate = map_lambda(sn.nodeparam{ist}{end}.setupTime);
267 alphascv = map_scv(sn.nodeparam{ist}{end}.setupTime);
268 betarate = map_lambda(sn.nodeparam{ist}{end}.delayoffTime);
269 betascv = map_scv(sn.nodeparam{ist}{end}.delayoffTime);
270 % Multiclass decomposition: compute per-class
271 % service rates and traffic intensities, solve
272 % aggregate setup-delayoff queue, then
273 % decompose per class proportionally
274 mu_k = zeros(1, K);
275 lambda_k = zeros(1, K);
276 active_k = false(1, K);
277 for k=1:K
278 pie_k = pie{ist}{k};
279 if ~isnan(pie_k(1))
280 mu_k(k) = 1 / (pie_k * inv(-D0{ist,k}) * ones(size(pie_k))');
281 c = find(sn.chains(:,k), 1);
282 lambda_k(k) = rates{ist,c}(k);
283 active_k(k) = true;
284 end
285 end
286 if any(active_k)
287 rho_k = lambda_k ./ mu_k;
288 rho_k(~active_k) = 0;
289 rho_total = sum(rho_k);
290 if rho_total > 0
291 % Aggregate service rate: weighted harmonic mean
292 aggrRate = aggrLambda / rho_total;
293 Q_total = qbd_setupdelayoff(aggrLambda, aggrRate, alpharate, alphascv, betarate, betascv);
294 for k=1:K
295 if active_k(k)
296 Qret{k} = Q_total * rho_k(k) / rho_total;
297 else
298 Qret{k} = NaN;
299 end
300 end
301 else
302 for k=1:K
303 Qret{k} = 0;
304 end
305 end
306 else
307 for k=1:K
308 Qret{k} = NaN;
309 end
310 end
311 else
312 [pdistr] = MMAPPH1FCFS(D, {pie{ist}{:}}, {D0{ist,:}}, 'ncDistr', maxLevel);
313 % rough approximation
314 for k=1:K
315 pdistr_k = abs(pdistr(1:(N(k)+1)));
316 pdistr_k(end) = abs(1-sum(pdistr(1:end-1)));
317 pdistr_k = pdistr_k / sum(pdistr_k(1:(N(k)+1)));
318 Qret{k} = max(0,min(N(k),(0:N(k))*pdistr_k(1:(N(k)+1))'));
319 end
320 end
321 end
322 end
323 else
324 for k=1:K
325 Qret{k} = sn.njobs(k);
326 end
327 end
328 end
329 QN(ist,:) = cell2mat(Qret);
330 for k=1:K
331 c = find(sn.chains(:,k),1);
332 TN(ist,k) = rates{ist,c}(k);
333 UN(ist,k) = TN(ist,k) * S(ist,k) / sn.nservers(ist);
334 QN(ist,k) = Qret{k};
335 if mapdcUsed
336 % For MAP/D/c, use exact results from qsys_mapdc
337 % No surrogate delay adjustment needed
338 RN(ist,k) = mapdcResult.meanSojournTime;
339 else
340 % add number of jobs at the surrogate delay server
341 QN(ist,k) = QN(ist,k) + TN(ist,k)*S(ist,k) * (sn.nservers(ist)-1)/sn.nservers(ist);
342 RN(ist,k) = QN(ist,k) ./ TN(ist,k);
343 end
344 end
345 end
346 end
347 else % not a station
348 switch sn.nodetype(ind)
349 case NodeType.Fork
350 % line_error(mfilename,'Fork nodes not supported yet by MAM solver.');
351 end
352 end
353 end
354 %it
355 %max(max(abs(TN-TN_1)))
356end
357totiter = it + 2;
358CN = sum(RN,1);
359QN = abs(QN);
360for it=1:2 % second pass to rescale again QN based on RN correction
361 for c=1:C
362 inchain = sn.inchain{c};
363 Nc = sum(sn.njobs(inchain));
364 if isfinite(Nc)
365 QNc = sum(sum(QN(:,inchain)));
366 QN(:,inchain) = QN(:,inchain) * (Nc / QNc);
367 end
368 for ind=1:I
369 for k=inchain
370 if sn.isstation(ind)
371 ist = sn.nodeToStation(ind);
372 % Skip stations using exact MAP/D/c solver (already have exact values)
373 if mapdcStations(ist)
374 continue;
375 end
376 if V(ist,k)>0
377 if isinf(sn.nservers(ist))
378 RN(ist,k) = S(ist,k);
379 else
380 RN(ist,k) = max([S(ist,k), QN(ist,k) ./ TN(ist,k)]);
381 end
382 else
383 RN(ist,k) = 0;
384 end
385 QN(ist,k) = RN(ist,k) .* TN(ist,k);
386 end
387 end
388 end
389 Nc = sum(sn.njobs(inchain)); % closed population
390 if Nc == 0 % if closed chain
391 QN(:,c)=0;
392 UN(:,c)=0;
393 RN(:,c)=0;
394 TN(:,c)=0;
395 CN(c)=0;
396 XN(c)=0;
397 end
398 end
399end
400end