LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
SolverENV.m
1classdef SolverENV < EnsembleSolver
2 % ENV - Ensemble environment solver for models with random environment changes
3 %
4 % SolverENV analyzes queueing networks operating in random environments where
5 % system parameters (arrival rates, service rates, routing) change according
6 % to an underlying environmental process. It solves ensemble models by analyzing
7 % each environmental stage and computing environment-averaged performance metrics.
8 %
9 % @brief Environment solver for networks in random changing environments
10 %
11 % Key characteristics:
12 % - Random environment with multiple operational stages
13 % - Environmental process governing parameter changes
14 % - Ensemble model analysis across environment stages
15 % - Environment-averaged performance computation
16 % - Stage-dependent system behavior modeling
17 %
18 % Environment solver features:
19 % - Multi-stage environmental modeling
20 % - Stage transition matrix analysis
21 % - Weighted performance metric computation
22 % - Environmental ensemble solution
23 % - Adaptive parameter modeling
24 %
25 % SolverENV is ideal for:
26 % - Systems with time-varying parameters
27 % - Networks subject to environmental fluctuations
28 % - Multi-mode operational system analysis
29 % - Performance under uncertainty modeling
30 % - Adaptive system behavior analysis
31 %
32 % Example:
33 % @code
34 % env_model = Environment(stages, transitions); % Define environment
35 % solver = SolverENV(env_model, @SolverMVA, options);
36 % metrics = solver.getEnsembleAvg(); % Environment-averaged metrics
37 % @endcode
38 %
39 % Copyright (c) 2012-2026, Imperial College London
40 % All rights reserved.
41
42
43 properties
44 env; % user-supplied representation of each stage transition
45 envObj;
46 sn;
47 resetFromMarginal;
48 resetEnvRates; % function implementing the reset policy for environment rates
49 stateDepMethod = ''; % state-dependent method configuration
50 SMPMethod = false; % Use DTMC-based computation for Semi-Markov Processes
51 % Enhanced init properties (aligned with JAR)
52 ServerNum; % Cell array of server counts per class
53 SRates; % Cell array of service rates per class
54 E0; % Rate matrix
55 Eutil; % Infinitesimal generator
56 transitionCdfs; % Transition CDF functions
57 sojournCdfs; % Sojourn time CDF functions
58 dtmcP; % DTMC transition matrix
59 holdTimeMatrix; % Hold time matrix
60 newMethod = false; % Use DTMC-based computation
61 compression = false; % Use Courtois decomposition
62 compressionResult; % Results from Courtois decomposition
63 Ecompress; % Number of macro-states after compression
64 MS; % Macro-state partition
65 end
66
67 methods
68 function self = SolverENV(renv, solverFactory, options)
69 % SELF = SOLVERENV(ENV,SOLVERFACTORY,OPTIONS)
70 self@EnsembleSolver(renv, mfilename);
71 if nargin>=3 %exist('options','var')
72 self.setOptions(options);
73 else
74 self.setOptions(SolverENV.defaultOptions);
75 end
76
77 % Enable SMP method if specified in options
78 if isfield(self.options, 'method') && strcmpi(self.options.method, 'smp')
79 self.SMPMethod = true;
80 line_debug('ENV solver: SMP method enabled via options.method=''smp''');
81 end
82
83 self.envObj = renv;
84 self.ensemble = renv.getEnsemble;
85 env = renv.getEnv;
86 self.env = env;
87 for e=1:length(self.env)
88 self.sn{e} = self.ensemble{e}.getStruct;
89 self.setSolver(solverFactory(self.ensemble{e}),e);
90 end
91
92 for e=1:length(self.env)
93 for h=1:length(self.env)
94 self.resetFromMarginal{e,h} = renv.resetFun{e,h};
95 end
96 end
97
98 for e=1:length(self.env)
99 for h=1:length(self.env)
100 self.resetEnvRates{e,h} = renv.resetEnvRatesFun{e,h};
101 end
102 end
103
104 for e=1:length(self.env)
105 for h=1:length(self.env)
106 if isa(self.env{e,h},'Disabled')
107 self.env{e,h} = Exp(0);
108 elseif ~isa(self.env{e,h},'Markovian') && ~self.SMPMethod
109 line_error(mfilename,sprintf('The distribution of the environment transition from stage %d to %d is not supported by the %s solver. Use method=''smp'' for non-Markovian distributions.',e,h,self.getName));
110 end
111 end
112 end
113
114 for e=1:length(self.ensemble)
115 if ~self.solvers{e}.supports(self.ensemble{e})
116 line_error(mfilename,sprintf('Model in the environment stage %d is not supported by the %s solver.',e,self.getName));
117 end
118 end
119 end
120
121 function setStateDepMethod(self, method)
122 % SETSTATEDEPMETHOD(METHOD) Sets the state-dependent method
123 if isempty(method)
124 line_error(mfilename, 'State-dependent method cannot be null or empty.');
125 end
126 self.stateDepMethod = method;
127 end
128
129 function setNewMethod(self, flag)
130 % SETNEWMETHOD(FLAG) Enable/disable DTMC-based computation
131 self.newMethod = flag;
132 end
133
134 function setCompression(self, flag)
135 % SETCOMPRESSION(FLAG) Enable/disable Courtois decomposition
136 self.compression = flag;
137 end
138
139 function [p, eps, epsMax, q] = ctmc_decompose(self, Q, MS)
140 % CTMC_DECOMPOSE Perform CTMC decomposition using configured method
141 % [p, eps, epsMax, q] = CTMC_DECOMPOSE(Q, MS)
142 %
143 % Uses options.config.decomp to select the decomposition algorithm:
144 % 'courtois' - Courtois decomposition (default)
145 % 'kms' - Koury-McAllister-Stewart method
146 % 'takahashi' - Takahashi's method
147 % 'multi' - Multigrid method (requires MSS)
148 %
149 % Returns:
150 % p - steady-state probability vector
151 % eps - NCD index
152 % epsMax - max acceptable eps value
153 % q - randomization coefficient
154
155 % Get decomposition/aggregation method from options
156 if isfield(self.options, 'config') && isfield(self.options.config, 'da')
157 method = self.options.config.da;
158 else
159 method = 'courtois';
160 end
161
162 % Get numsteps for iterative methods
163 if isfield(self.options, 'config') && isfield(self.options.config, 'da_iter')
164 numsteps = self.options.config.da_iter;
165 else
166 numsteps = 10;
167 end
168
169 switch lower(method)
170 case 'courtois'
171 [p, ~, ~, eps, epsMax, ~, ~, ~, q] = ctmc_courtois(Q, MS);
172 case 'kms'
173 [p, ~, ~, eps, epsMax] = ctmc_kms(Q, MS, numsteps);
174 q = 1.05 * max(max(abs(Q)));
175 case 'takahashi'
176 [p, ~, ~, ~, eps, epsMax] = ctmc_takahashi(Q, MS, numsteps);
177 q = 1.05 * max(max(abs(Q)));
178 case 'multi'
179 % Multi requires MSS (macro-macro-states), default to singletons
180 nMacro = size(MS, 1);
181 MSS = cell(nMacro, 1);
182 for i = 1:nMacro
183 MSS{i} = i;
184 end
185 [p, ~, ~, ~, eps, epsMax] = ctmc_multi(Q, MS, MSS);
186 q = 1.05 * max(max(abs(Q)));
187 otherwise
188 line_error(mfilename, sprintf('Unknown decomposition method: %s', method));
189 end
190 end
191
192 function bool = converged(self, it) % convergence test at iteration it
193 % BOOL = CONVERGED(IT) % CONVERGENCE TEST AT ITERATION IT
194 % Computes max relative absolute difference of queue lengths between iterations
195 % Aligned with JAR SolverEnv.converged() implementation
196
197 bool = false;
198 if it <= 1
199 return
200 end
201
202 E = self.getNumberOfModels;
203 M = self.ensemble{1}.getNumberOfStations;
204 K = self.ensemble{1}.getNumberOfClasses;
205
206 % Check convergence per class (aligned with JAR structure)
207 for k = 1:K
208 % Build QEntry and QExit matrices (M x E) for this class
209 QEntry = zeros(M, E);
210 QExit = zeros(M, E);
211 for e = 1:E
212 % Skip stages where analysis failed
213 res_curr = self.results{it,e};
214 res_prev = self.results{it-1,e};
215 for i = 1:M
216 if ~isempty(res_curr) && isfield(res_curr, 'Tran') ...
217 && isstruct(res_curr.Tran.Avg) && isfield(res_curr.Tran.Avg, 'Q')
218 Qik_curr = res_curr.Tran.Avg.Q{i,k};
219 if isstruct(Qik_curr) && isfield(Qik_curr, 'metric') && ~isempty(Qik_curr.metric)
220 QExit(i,e) = Qik_curr.metric(1);
221 end
222 end
223 if ~isempty(res_prev) && isfield(res_prev, 'Tran') ...
224 && isstruct(res_prev.Tran.Avg) && isfield(res_prev.Tran.Avg, 'Q')
225 Qik_prev = res_prev.Tran.Avg.Q{i,k};
226 if isstruct(Qik_prev) && isfield(Qik_prev, 'metric') && ~isempty(Qik_prev.metric)
227 QEntry(i,e) = Qik_prev.metric(1);
228 end
229 end
230 end
231 end
232
233 % Compute max relative absolute difference using maxpe
234 % maxpe computes max(abs(1 - approx./exact)) = max(abs((approx-exact)./exact))
235 % This matches JAR's Matrix.maxAbsDiff() implementation
236 maxDiff = maxpe(QExit(:), QEntry(:));
237 if isempty(maxDiff)
238 maxDiff = 0; % Handle case where all QEntry values are zero
239 end
240 if isnan(maxDiff) || isinf(maxDiff)
241 return % Non-convergence on invalid values
242 end
243 if maxDiff >= self.options.iter_tol
244 return % Not converged
245 end
246 end
247 bool = true;
248 end
249
250
251 function runAnalyzer(self)
252 % RUNANALYZER()
253 % Run the ensemble solver iteration
254 line_debug('ENV solver starting: nstages=%d, method=%s', self.getNumberOfModels, self.options.method);
255
256 % Show library attribution if verbose and not yet shown
257 if self.options.verbose ~= VerboseLevel.SILENT && ~GlobalConstants.isLibraryAttributionShown()
258 libs = SolverENV.getLibrariesUsed([], self.options);
259 if ~isempty(libs)
260 line_printf('The solver will leverage %s.\n', strjoin(libs, ', '));
261 GlobalConstants.setLibraryAttributionShown(true);
262 end
263 end
264
265 iterate(self);
266 end
267
268 function init(self)
269 % INIT()
270 % Initialize the environment solver with enhanced data structures
271 % aligned with JAR SolverEnv implementation
272 line_debug('ENV solver init: initializing environment data structures');
273 options = self.options;
274 if isfield(options,'seed')
275 Solver.resetRandomGeneratorSeed(options.seed);
276 end
277 self.envObj.init();
278
279 % Initialize ServerNum and SRates matrices (aligned with JAR)
280 E = self.getNumberOfModels;
281 M = self.sn{1}.nstations;
282 K = self.sn{1}.nclasses;
283
284 self.ServerNum = cell(K, 1);
285 self.SRates = cell(K, 1);
286 for k = 1:K
287 self.ServerNum{k} = zeros(M, E);
288 self.SRates{k} = zeros(M, E);
289 for e = 1:E
290 for m = 1:M
291 self.ServerNum{k}(m, e) = self.sn{e}.nservers(m);
292 self.SRates{k}(m, e) = self.sn{e}.rates(m, k);
293 end
294 end
295 end
296
297 % Build rate matrix E0
298 self.E0 = zeros(E, E);
299 for e = 1:E
300 for h = 1:E
301 if ~isa(self.envObj.env{e,h}, 'Disabled')
302 self.E0(e, h) = self.envObj.env{e,h}.getRate();
303 end
304 end
305 end
306 self.Eutil = ctmc_makeinfgen(self.E0);
307
308 % Initialize transition CDFs
309 self.transitionCdfs = cell(E, E);
310 for e = 1:E
311 for h = 1:E
312 if ~isa(self.envObj.env{e,h}, 'Disabled')
313 envDist = self.envObj.env{e,h};
314 self.transitionCdfs{e,h} = @(t) envDist.evalCDF(t);
315 else
316 self.transitionCdfs{e,h} = @(t) 0;
317 end
318 end
319 end
320
321 % Initialize sojourn CDFs
322 self.sojournCdfs = cell(E, 1);
323 for e = 1:E
324 self.sojournCdfs{e} = @(t) self.computeSojournCdf(e, t);
325 end
326
327 % newMethod: Use DTMC-based computation instead of CTMC
328 % Verified numerical integration for Semi-Markov Process DTMC transition probabilities
329 if self.newMethod
330 line_debug('ENV using DTMC-based computation (newMethod=true)');
331 self.dtmcP = zeros(E, E);
332 for k = 1:E
333 for e = 1:E
334 if k == e || isa(self.envObj.env{k,e}, 'Disabled')
335 self.dtmcP(k, e) = 0.0;
336 else
337 % Compute the upper limit of the sojourn time
338 epsilon = 1e-8;
339 T = 1;
340 while self.transitionCdfs{k,e}(T) < 1.0 - epsilon
341 T = T * 2;
342 if T > 1e6 % safety limit
343 break;
344 end
345 end
346 % Adaptive number of integration intervals based on T
347 N = max(1000, round(T * 100));
348 dt = T / N;
349 sumVal = 0;
350 for i = 0:(N-1)
351 t0 = i * dt;
352 t1 = t0 + dt;
353 deltaF = self.transitionCdfs{k,e}(t1) - self.transitionCdfs{k,e}(t0);
354 survival = 1;
355 for h = 1:E
356 if h ~= k && h ~= e && ~isa(self.envObj.env{k,h}, 'Disabled')
357 % Use midpoint for better accuracy
358 tmid = (t0 + t1) / 2.0;
359 survival = survival * (1.0 - self.envObj.env{k,h}.evalCDF(tmid));
360 end
361 end
362 sumVal = sumVal + deltaF * survival;
363 end
364 self.dtmcP(k, e) = sumVal;
365 end
366 end
367 end
368
369 % Solve DTMC for stationary distribution
370 dtmcPie = dtmc_solve(self.dtmcP);
371
372 % Calculate hold times using numerical integration
373 self.holdTimeMatrix = self.computeHoldTime(E);
374
375 % Compute steady-state probabilities
376 pi = zeros(1, E);
377 denomSum = 0;
378 for e = 1:E
379 denomSum = denomSum + dtmcPie(e) * self.holdTimeMatrix(e);
380 end
381 for k = 1:E
382 pi(k) = dtmcPie(k) * self.holdTimeMatrix(k) / denomSum;
383 end
384 self.envObj.probEnv = pi;
385
386 % Update embedding weights
387 newEmbweight = zeros(E, E);
388 for e = 1:E
389 sumVal = 0.0;
390 for h = 1:E
391 if h ~= e
392 sumVal = sumVal + pi(h) * self.E0(h, e);
393 end
394 end
395 for k = 1:E
396 if k == e
397 newEmbweight(k, e) = 0;
398 else
399 if sumVal > 0
400 newEmbweight(k, e) = pi(k) * self.E0(k, e) / sumVal;
401 end
402 end
403 end
404 end
405 self.envObj.probOrig = newEmbweight;
406 end
407
408 % Compression: Use Courtois decomposition for large environments
409 if self.compression
410 line_debug('ENV using compression (Courtois decomposition)');
411 self.applyCompression(E, M, K);
412 end
413 end
414
415 function applyCompression(self, E, M, K)
416 % APPLYCOMPRESSION Apply Courtois decomposition to reduce environment size
417 % This method finds a good partition of the environment states and
418 % creates compressed macro-state networks.
419
420 % Find best partition
421 if E <= 10
422 self.MS = self.findBestPartition(E);
423 else
424 % Beam search for large environments
425 self.MS = self.beamSearchPartition(E);
426 end
427
428 if isempty(self.MS)
429 % No compression possible, use singletons
430 self.MS = cell(E, 1);
431 for i = 1:E
432 self.MS{i} = i;
433 end
434 self.Ecompress = E;
435 return;
436 end
437
438 self.Ecompress = length(self.MS);
439
440 % Apply decomposition/aggregation
441 [p, eps, epsMax, q] = self.ctmc_decompose(self.Eutil, self.MS);
442
443 if eps > epsMax
444 line_warning(mfilename, 'Environment cannot be effectively compressed (eps > epsMax).');
445 end
446
447 % Store compression results
448 self.compressionResult.p = p;
449 self.compressionResult.eps = eps;
450 self.compressionResult.epsMax = epsMax;
451 self.compressionResult.q = q;
452
453 % Update probEnv with macro-state probabilities
454 pMacro = zeros(1, self.Ecompress);
455 for i = 1:self.Ecompress
456 pMacro(i) = sum(p(self.MS{i}));
457 end
458 self.envObj.probEnv = pMacro;
459
460 % Compute micro-state probabilities within each macro-state
461 pmicro = zeros(E, 1);
462 for i = 1:self.Ecompress
463 blockProb = p(self.MS{i});
464 if sum(blockProb) > 0
465 pmicro(self.MS{i}) = blockProb / sum(blockProb);
466 end
467 end
468 self.compressionResult.pmicro = pmicro;
469 self.compressionResult.pMacro = pMacro;
470
471 % Update embedding weights for macro-states
472 Ecomp = self.Ecompress;
473 newEmbweight = zeros(Ecomp, Ecomp);
474 for e = 1:Ecomp
475 sumVal = 0.0;
476 for h = 1:Ecomp
477 if h ~= e
478 sumVal = sumVal + pMacro(h) * self.computeMacroRate(h, e);
479 end
480 end
481 for k = 1:Ecomp
482 if k == e
483 newEmbweight(k, e) = 0;
484 elseif sumVal > 0
485 newEmbweight(k, e) = pMacro(k) * self.computeMacroRate(k, e) / sumVal;
486 end
487 end
488 end
489 self.envObj.probOrig = newEmbweight;
490
491 % Build macro-state networks with weighted-average rates
492 macroEnsemble = cell(self.Ecompress, 1);
493 macroSolvers = cell(self.Ecompress, 1);
494 macroSn = cell(self.Ecompress, 1);
495
496 for i = 1:self.Ecompress
497 % Copy the first micro-state network
498 firstMicro = self.MS{i}(1);
499 macroEnsemble{i} = self.ensemble{firstMicro}.copy();
500
501 % Compute weighted-average rates
502 for m = 1:M
503 for k = 1:K
504 rateSum = 0;
505 for r = 1:length(self.MS{i})
506 microIdx = self.MS{i}(r);
507 w = pmicro(microIdx);
508 rateSum = rateSum + w * self.sn{microIdx}.rates(m, k);
509 end
510
511 % Update service rate
512 jobclass = macroEnsemble{i}.classes{k};
513 station = macroEnsemble{i}.stations{m};
514 if isa(station, 'Queue') || isa(station, 'Delay')
515 if rateSum > 0
516 station.setService(jobclass, Exp(rateSum));
517 end
518 end
519 end
520 end
521
522 macroEnsemble{i}.refreshStruct(true);
523 macroSn{i} = macroEnsemble{i}.getStruct(true);
524
525 % Create solver for macro-state
526 % Use the same solver factory pattern as original
527 macroSolvers{i} = SolverFluid(macroEnsemble{i}, self.solvers{firstMicro}.options);
528 end
529
530 % Replace ensemble and solvers with compressed versions
531 self.ensemble = macroEnsemble;
532 self.solvers = macroSolvers;
533 self.sn = macroSn;
534 end
535
536 function rate = computeMacroRate(self, fromMacro, toMacro)
537 % COMPUTEMACRORATE Compute transition rate between macro-states
538 rate = 0;
539 for i = 1:length(self.MS{fromMacro})
540 mi = self.MS{fromMacro}(i);
541 for j = 1:length(self.MS{toMacro})
542 mj = self.MS{toMacro}(j);
543 rate = rate + self.compressionResult.pmicro(mi) * self.E0(mi, mj);
544 end
545 end
546 end
547
548 function MS = findBestPartition(self, E)
549 % FINDBESTPARTITION Find the best partition for small environments (E <= 10)
550 % Uses exhaustive search over all possible partitions
551
552 % Start with singletons
553 bestMS = cell(E, 1);
554 for i = 1:E
555 bestMS{i} = i;
556 end
557
558 [~, bestEps, bestEpsMax] = self.ctmc_decompose(self.Eutil, bestMS);
559 if isempty(bestEps) || isnan(bestEps)
560 MS = bestMS;
561 return;
562 end
563
564 % Try merging pairs
565 for i = 1:E
566 for j = (i+1):E
567 testMS = cell(E-1, 1);
568 idx = 1;
569 for k = 1:E
570 if k == i
571 testMS{idx} = [i, j];
572 idx = idx + 1;
573 elseif k ~= j
574 testMS{idx} = k;
575 idx = idx + 1;
576 end
577 end
578
579 [~, testEps, testEpsMax] = self.ctmc_decompose(self.Eutil, testMS);
580 if ~isempty(testEps) && ~isnan(testEps) && testEps < bestEps
581 bestEps = testEps;
582 bestEpsMax = testEpsMax;
583 bestMS = testMS;
584 end
585 end
586 end
587
588 MS = bestMS;
589 self.Ecompress = length(MS);
590 end
591
592 function MS = beamSearchPartition(self, E)
593 % BEAMSEARCHPARTITION Beam search for large environments (E > 10)
594
595 B = 3; % Beam width
596 alpha = 0.01; % Coupling threshold
597
598 if isfield(self.options, 'config') && isfield(self.options.config, 'env_alpha')
599 alpha = self.options.config.env_alpha;
600 end
601
602 % Initialize with singletons
603 singletons = cell(E, 1);
604 for i = 1:E
605 singletons{i} = i;
606 end
607 beam = {singletons};
608
609 bestSeen = singletons;
610 [~, bestEps] = self.ctmc_decompose(self.Eutil, bestSeen);
611
612 % Iteratively merge blocks
613 for depth = 1:(E-1)
614 candidates = {};
615
616 for b = 1:length(beam)
617 ms = beam{b};
618 nBlocks = length(ms);
619
620 % Try all pairwise merges
621 for i = 1:nBlocks
622 for j = (i+1):nBlocks
623 % Create merged partition
624 trial = cell(nBlocks - 1, 1);
625 idx = 1;
626 for k = 1:nBlocks
627 if k == i
628 trial{idx} = [ms{i}(:); ms{j}(:)];
629 idx = idx + 1;
630 elseif k ~= j
631 trial{idx} = ms{k};
632 idx = idx + 1;
633 end
634 end
635
636 [~, childEps, childEpsMax] = self.ctmc_decompose(self.Eutil, trial);
637
638 if ~isempty(childEps) && ~isnan(childEps) && childEps > 0
639 cost = childEps - childEpsMax + alpha * depth;
640 candidates{end+1} = {trial, cost};
641
642 if cost < bestEps
643 bestEps = cost;
644 bestSeen = trial;
645 end
646 end
647 end
648 end
649 end
650
651 if isempty(candidates)
652 break;
653 end
654
655 % Sort by cost and keep top B
656 costs = cellfun(@(x) x{2}, candidates);
657 [~, sortIdx] = sort(costs);
658 beam = {};
659 for i = 1:min(B, length(sortIdx))
660 beam{end+1} = candidates{sortIdx(i)}{1};
661 end
662 end
663
664 MS = bestSeen;
665 self.Ecompress = length(MS);
666 end
667
668 function holdTime = computeHoldTime(self, E)
669 % COMPUTEHOLDTIME Compute expected holding times using numerical integration
670 holdTime = zeros(1, E);
671 for k = 1:E
672 % Survival function: 1 - sojournCDF
673 surv = @(t) 1 - self.sojournCdfs{k}(t);
674
675 % Compute upper limit
676 upperLimit = 10;
677 while surv(upperLimit) > 1e-8
678 upperLimit = upperLimit * 2;
679 if upperLimit > 1e6 % safety limit
680 break;
681 end
682 end
683
684 % Simpson's rule integration
685 N = 10000;
686 dt = upperLimit / N;
687 integral = 0;
688 for i = 0:(N-1)
689 t0 = i * dt;
690 t1 = t0 + dt;
691 tmid = (t0 + t1) / 2;
692 % Simpson's rule: (f(a) + 4*f(mid) + f(b)) * h/6
693 integral = integral + (surv(t0) + 4*surv(tmid) + surv(t1)) * dt / 6;
694 end
695 holdTime(k) = integral;
696 end
697 end
698
699 function cdf = computeSojournCdf(self, e, t)
700 % COMPUTESOJOURNCDF Compute sojourn time CDF for environment stage e
701 E = self.getNumberOfModels;
702 surv = 1.0;
703 for h = 1:E
704 if h ~= e
705 surv = surv * (1 - self.transitionCdfs{e,h}(t));
706 end
707 end
708 cdf = 1 - surv;
709 end
710
711 function pre(self, it)
712 % PRE(IT)
713
714 if it==1
715 for e=list(self)
716 try
717 if isinf(self.getSolver(e).options.timespan(2))
718 [QN,~,~,~] = self.getSolver(e).getAvg();
719 else
720 [QNt,~,~] = self.getSolver(e).getTranAvg();
721 % Handle case where getTranAvg returns NaN for disabled/empty results
722 QN = zeros(size(QNt));
723 for i = 1:size(QNt,1)
724 for k = 1:size(QNt,2)
725 if isstruct(QNt{i,k}) && isfield(QNt{i,k}, 'metric')
726 QN(i,k) = QNt{i,k}.metric(end);
727 else
728 QN(i,k) = 0; % Default to 0 for NaN/invalid entries
729 end
730 end
731 end
732 end
733 if ~isa(self.solvers{e}, 'SolverFluid')
734 QN = self.roundMarginalForDiscreteSolver(QN, self.sn{e});
735 end
736 self.ensemble{e}.initFromMarginal(QN);
737 catch
738 % Skip pre-initialization for stages where getAvg fails
739 end
740 end
741 end
742 end
743
744 % solves model in stage e
745 function [results_e, runtime] = analyze(self, it, e)
746 % [RESULTS_E, RUNTIME] = ANALYZE(IT, E)
747 results_e = struct();
748 results_e.('Tran') = struct();
749 results_e.Tran.('Avg') = [];
750 T0 = tic;
751 runtime = toc(T0);
752 %% initialize
753 try
754 [Qt,Ut,Tt] = self.ensemble{e}.getTranHandles;
755 self.solvers{e}.reset();
756 [QNt,UNt,TNt] = self.solvers{e}.getTranAvg(Qt,Ut,Tt);
757 results_e.Tran.Avg.Q = QNt;
758 results_e.Tran.Avg.U = UNt;
759 results_e.Tran.Avg.T = TNt;
760 catch
761 % Skip analysis for stages where transient analysis fails
762 end
763 end
764
765
766 function post(self, it)
767 % POST(IT)
768
769 E = self.getNumberOfModels;
770 M = self.ensemble{1}.getNumberOfStations;
771 K = self.ensemble{1}.getNumberOfClasses;
772 for e=1:E
773 if isempty(self.results{it,e}) || ~isfield(self.results{it,e}, 'Tran') ...
774 || ~isstruct(self.results{it,e}.Tran.Avg) || ~isfield(self.results{it,e}.Tran.Avg, 'Q')
775 for h = 1:E
776 Qexit{e,h} = zeros(M, K);
777 Uexit{e,h} = zeros(M, K);
778 Texit{e,h} = zeros(M, K);
779 end
780 continue
781 end
782 for h = 1:E
783 Qexit{e,h} = zeros(size(self.results{it,e}.Tran.Avg.Q));
784 Uexit{e,h} = zeros(size(self.results{it,e}.Tran.Avg.U));
785 Texit{e,h} = zeros(size(self.results{it,e}.Tran.Avg.T));
786 for i=1:size(self.results{it,e}.Tran.Avg.Q,1)
787 for r=1:size(self.results{it,e}.Tran.Avg.Q,2)
788 Qir = self.results{it,e}.Tran.Avg.Q{i,r};
789 % Check if result is a valid struct with required fields
790 if isstruct(Qir) && isfield(Qir, 't') && isfield(Qir, 'metric') && ~isempty(Qir.t)
791 w{e,h} = [0, map_cdf(self.envObj.proc{e}{h}, Qir.t(2:end)) - map_cdf(self.envObj.proc{e}{h}, Qir.t(1:end-1))]';
792 if ~isnan(w{e,h})
793 Qexit{e,h}(i,r) = Qir.metric'*w{e,h}/sum(w{e,h});
794 Uir = self.results{it,e}.Tran.Avg.U{i,r};
795 if isstruct(Uir) && isfield(Uir, 'metric') && ~isempty(Uir.metric)
796 Uexit{e,h}(i,r) = Uir.metric'*w{e,h}/sum(w{e,h});
797 end
798 Tir = self.results{it,e}.Tran.Avg.T{i,r};
799 if isstruct(Tir) && isfield(Tir, 'metric') && ~isempty(Tir.metric)
800 Texit{e,h}(i,r) = Tir.metric'*w{e,h}/sum(w{e,h});
801 end
802 else
803 Qexit{e,h}(i,r) = 0;
804 Uexit{e,h}(i,r) = 0;
805 Texit{e,h}(i,r) = 0;
806 end
807 else
808 w{e,h} = 0;
809 end
810 end
811 end
812 end
813 end
814
815 Qentry = cell(1,E); % average entry queue-length
816 for e = 1:E
817 % Skip stages where analysis failed (no valid Tran results)
818 if isempty(self.results{it,e}) || ~isfield(self.results{it,e}, 'Tran') ...
819 || ~isstruct(self.results{it,e}.Tran.Avg) || ~isfield(self.results{it,e}.Tran.Avg, 'Q')
820 continue
821 end
822 Qentry{e} = zeros(size(Qexit{e}));
823 for h=1:E
824 % probability of coming from h to e \times resetFun(Qexit from h to e
825 if self.envObj.probOrig(h,e) > 0
826 Qentry{e} = Qentry{e} + self.envObj.probOrig(h,e) * self.resetFromMarginal{h,e}(Qexit{h,e});
827 end
828 end
829 if ~isa(self.solvers{e}, 'SolverFluid')
830 Qentry{e} = self.roundMarginalForDiscreteSolver(Qentry{e}, self.sn{e});
831 end
832 self.solvers{e}.reset();
833 self.ensemble{e}.initFromMarginal(Qentry{e});
834 end
835
836 % Update transition rates between stages if state-dependent method
837 if isfield(self.options, 'method') && strcmp(self.options.method, 'statedep')
838 for e = 1:E
839 for h = 1:E
840 if ~isa(self.envObj.env{e,h}, 'Disabled') && ~isempty(self.resetEnvRates{e,h})
841 self.envObj.env{e,h} = self.resetEnvRates{e,h}(...
842 self.envObj.env{e,h}, Qexit{e,h}, Uexit{e,h}, Texit{e,h});
843 end
844 end
845 end
846 % Reinitialize environment after rate updates
847 self.envObj.init();
848 end
849 end
850
851 function finish(self)
852 % FINISH()
853
854 it = size(self.results,1); % use last iteration
855 E = self.getNumberOfModels;
856 M = self.ensemble{1}.getNumberOfStations;
857 K = self.ensemble{1}.getNumberOfClasses;
858 for e=1:E
859 QExit{e}=zeros(M, K);
860 UExit{e}=zeros(M, K);
861 TExit{e}=zeros(M, K);
862 if it>0 && ~isempty(self.results{it,e}) && isfield(self.results{it,e}, 'Tran') ...
863 && isstruct(self.results{it,e}.Tran.Avg) && isfield(self.results{it,e}.Tran.Avg, 'Q')
864 for i=1:size(self.results{it,e}.Tran.Avg.Q,1)
865 for r=1:size(self.results{it,e}.Tran.Avg.Q,2)
866 Qir = self.results{it,e}.Tran.Avg.Q{i,r};
867 % Check if result is a valid struct with required fields
868 if isstruct(Qir) && isfield(Qir, 't') && isfield(Qir, 'metric') && ~isempty(Qir.t)
869 w{e} = [0, map_cdf(self.envObj.holdTime{e}, Qir.t(2:end)) - map_cdf(self.envObj.holdTime{e}, Qir.t(1:end-1))]';
870 QExit{e}(i,r) = Qir.metric'*w{e}/sum(w{e});
871 Uir = self.results{it,e}.Tran.Avg.U{i,r};
872 if isstruct(Uir) && isfield(Uir, 'metric') && ~isempty(Uir.metric)
873 UExit{e}(i,r) = Uir.metric'*w{e}/sum(w{e});
874 else
875 UExit{e}(i,r) = 0;
876 end
877 Tir = self.results{it,e}.Tran.Avg.T{i,r};
878 if isstruct(Tir) && isfield(Tir, 'metric') && ~isempty(Tir.metric)
879 TExit{e}(i,r) = Tir.metric'*w{e}/sum(w{e});
880 else
881 TExit{e}(i,r) = 0;
882 end
883 else
884 QExit{e}(i,r) = 0;
885 UExit{e}(i,r) = 0;
886 TExit{e}(i,r) = 0;
887 end
888 end
889 end
890 end
891 % for h = 1:E
892 % QE{e,h} = zeros(size(self.results{it,e}.Tran.Avg.Q));
893 % for i=1:size(self.results{it,e}.Tran.Avg.Q,1)
894 % for r=1:size(self.results{it,e}.Tran.Avg.Q,2)
895 % w{e,h} = [0, map_cdf(self.envObj.proc{e}{h}, self.results{it,e}.Tran.Avg.Q{i,r}(2:end,2)) - map_cdf(self.envObj.proc{e}{h}, self.results{it,e}.Tran.Avg.Q{i,r}(1:end-1,2))]';
896 % if ~isnan(w{e,h})
897 % QE{e,h}(i,r) = self.results{it,e}.Tran.Avg.Q{i,r}(:,1)'*w{e,h}/sum(w{e,h});
898 % else
899 % QE{e,h}(i,r) = 0;
900 % end
901 % end
902 % end
903 % end
904 end
905
906 Qval=0*QExit{e};
907 Uval=0*UExit{e};
908 Tval=0*TExit{e};
909 for e=1:E
910 Qval = Qval + self.envObj.probEnv(e) * QExit{e}; % to check
911 Uval = Uval + self.envObj.probEnv(e) * UExit{e}; % to check
912 Tval = Tval + self.envObj.probEnv(e) * TExit{e}; % to check
913 end
914 self.result.Avg.Q = Qval;
915 % self.result.Avg.R = R;
916 % self.result.Avg.X = X;
917 self.result.Avg.U = Uval;
918 self.result.Avg.T = Tval;
919 % self.result.Avg.C = C;
920 %self.result.runtime = runtime;
921 %if self.options.verbose
922 % line_printf('\n');
923 %end
924 end
925
926 function name = getName(self)
927 % NAME = GETNAME()
928
929 name = mfilename;
930 end
931
932 function [renvInfGen, stageInfGen, renvEventFilt, stageEventFilt, renvEvents, stageEvents] = getGenerator(self)
933 % [renvInfGen, stageInfGen, renvEventFilt, stageEventFilt, renvEvents, stageEvents] = getGenerator(self)
934 %
935 % Returns the infinitesimal generator matrices for the random environment model.
936 %
937 % Outputs:
938 % renvInfGen - Combined infinitesimal generator for the random environment (flattened)
939 % stageInfGen - Cell array of infinitesimal generators for each stage
940 % renvEventFilt - Cell array (E x E) of event filtration matrices for environment transitions
941 % stageEventFilt - Cell array of event filtration matrices for each stage
942 % renvEvents - Cell array of Event objects for environment transitions
943 % stageEvents - Cell array of synchronization maps for each stage
944
945 E = self.getNumberOfModels;
946 stageInfGen = cell(1,E);
947 stageEventFilt = cell(1,E);
948 stageEvents = cell(1,E);
949 for e=1:E
950 if isa(self.solvers{e},'SolverCTMC')
951 [stageInfGen{e}, stageEventFilt{e}, stageEvents{e}] = self.solvers{e}.getGenerator();
952 else
953 line_error(mfilename,'This method requires SolverENV to be instantiated with the CTMC solver.');
954 end
955 end
956
957 % Get number of states for each stage
958 nstates = cellfun(@(g) size(g, 1), stageInfGen);
959
960 % Get number of phases for each transition distribution
961 nphases = zeros(E, E);
962 for i=1:E
963 for j=1:E
964 if ~isempty(self.env{i,j}) && ~isa(self.env{i,j}, 'Disabled')
965 nphases(i,j) = self.env{i,j}.getNumberOfPhases();
966 else
967 nphases(i,j) = 1;
968 end
969 end
970 end
971 % Adjust diagonal (self-transitions have one less phase in the Kronecker expansion)
972 nphases = nphases - eye(E);
973
974 % Initialize block cell structure for the random environment generator
975 renvInfGen = cell(E,E);
976
977 for e=1:E
978 % Diagonal block: stage infinitesimal generator
979 renvInfGen{e,e} = stageInfGen{e};
980 for h=1:E
981 if h~=e
982 % Off-diagonal blocks: reset matrices (identity with appropriate dimensions)
983 minStates = min(nstates(e), nstates(h));
984 resetMatrix_eh = sparse(nstates(e), nstates(h));
985 for i=1:minStates
986 resetMatrix_eh(i,i) = 1;
987 end
988 renvInfGen{e,h} = resetMatrix_eh;
989 end
990 end
991 end
992
993 % Build environment transition events and expand generator with phase structure
994 renvEvents = cell(1,0);
995 for e=1:E
996 for h=1:E
997 if h~=e
998 % Get D0 (phase generator) for transition from e to h
999 if isempty(self.env{e,h}) || isa(self.env{e,h}, 'Disabled')
1000 D0 = zeros(0);
1001 else
1002 proc = self.env{e,h}.getProcess();
1003 D0 = proc{1};
1004 end
1005 % Kronecker sum with diagonal block
1006 renvInfGen{e,e} = krons(renvInfGen{e,e}, D0);
1007
1008 % Get D1 (completion rate matrix) and initial probability vector pie
1009 if isempty(self.env{h,e}) || isa(self.env{h,e}, 'Disabled') || any(isnan(map_pie(self.env{h,e}.getProcess())))
1010 pie = ones(1, nphases(h,e));
1011 else
1012 pie = map_pie(self.env{h,e}.getProcess());
1013 end
1014
1015 if isempty(self.env{e,h}) || isa(self.env{e,h}, 'Disabled')
1016 D1 = zeros(0);
1017 else
1018 proc = self.env{e,h}.getProcess();
1019 D1 = proc{2};
1020 end
1021
1022 % Kronecker product for off-diagonal block
1023 onePhase = ones(nphases(e,h), 1);
1024 kronArg = D1 * onePhase * pie;
1025 renvInfGen{e,h} = kron(renvInfGen{e,h}, sparse(kronArg));
1026
1027 % Create environment transition events
1028 for i=1:self.ensemble{e}.getNumberOfNodes
1029 renvEvents{1,end+1} = Event(EventType.STAGE, i, NaN, NaN, [e,h]); %#ok<AGROW>
1030 end
1031
1032 % Handle other stages (f != e, f != h)
1033 for f=1:E
1034 if f~=h && f~=e
1035 if isempty(self.env{f,h}) || isa(self.env{f,h}, 'Disabled') || any(isnan(map_pie(self.env{f,h}.getProcess())))
1036 pie_fh = ones(1, nphases(f,h));
1037 else
1038 pie_fh = map_pie(self.env{f,h}.getProcess());
1039 end
1040 oneVec = ones(nphases(e,h), 1);
1041 renvInfGen{e,f} = kron(renvInfGen{e,f}, oneVec * pie_fh);
1042 end
1043 end
1044 end
1045 end
1046 end
1047
1048 % Build event filtration matrices for environment transitions
1049 % Each renvEventFilt{e,h} isolates transitions from stage e to stage h
1050 renvEventFilt = cell(E,E);
1051 for e=1:E
1052 for h=1:E
1053 tmpCell = cell(E,E);
1054 for e1=1:E
1055 for h1=1:E
1056 tmpCell{e1,h1} = renvInfGen{e1,h1};
1057 end
1058 end
1059 % Zero out diagonal blocks (internal stage transitions)
1060 for e1=1:E
1061 tmpCell{e1,e1} = tmpCell{e1,e1} * 0;
1062 end
1063 % Zero out off-diagonal blocks that don't match (e,h) transition
1064 for e1=1:E
1065 for h1=1:E
1066 if e~=e1 || h~=h1
1067 if e1~=h1 % Only zero out off-diagonal entries
1068 tmpCell{e1,h1} = tmpCell{e1,h1} * 0;
1069 end
1070 end
1071 end
1072 end
1073 renvEventFilt{e,h} = cell2mat(tmpCell);
1074 end
1075 end
1076
1077 % Flatten block structure into single matrix and normalize
1078 renvInfGen = cell2mat(renvInfGen);
1079 renvInfGen = ctmc_makeinfgen(renvInfGen);
1080 end
1081
1082 function varargout = getAvg(varargin)
1083 % [QNCLASS, UNCLASS, TNCLASS] = GETAVG()
1084 [varargout{1:nargout}] = getEnsembleAvg( varargin{:} );
1085 end
1086
1087 function [QNclass, UNclass, RNclass, TNclass, ANclass, WNclass] = getEnsembleAvg(self)
1088 % [QNCLASS, UNCLASS, TNCLASS] = GETENSEMBLEAVG()
1089 RNclass=[];
1090 ANclass=[];
1091 WNclass=[];
1092
1093 if isempty(self.result) || (isfield(self.options,'force') && self.options.force)
1094 iterate(self);
1095 if isempty(self.result)
1096 QNclass=[];
1097 UNclass=[];
1098 TNclass=[];
1099 return
1100 end
1101 end
1102 QNclass = self.result.Avg.Q;
1103 UNclass = self.result.Avg.U;
1104 TNclass = self.result.Avg.T;
1105 WNclass = QNclass ./ TNclass;
1106 RNclass = NaN*WNclass;
1107 ANclass = NaN*TNclass;
1108 end
1109
1110 function [AvgTable,QT,UT,TT] = getAvgTable(self,keepDisabled)
1111 % [AVGTABLE,QT,UT,TT] = GETAVGTABLE(SELF,KEEPDISABLED)
1112 % Return table of average station metrics
1113
1114 if nargin<2 %if ~exist('keepDisabled','var')
1115 keepDisabled = false;
1116 end
1117
1118 [QN,UN,~,TN] = getAvg(self);
1119 M = size(QN,1);
1120 K = size(QN,2);
1121 Q = self.result.Avg.Q;
1122 U = self.result.Avg.U;
1123 T = self.result.Avg.T;
1124 if isempty(QN)
1125 AvgTable = Table();
1126 QT = Table();
1127 UT = Table();
1128 TT = Table();
1129 elseif ~keepDisabled
1130 Qval = []; Uval = []; Tval = [];
1131 JobClass = {};
1132 Station = {};
1133 for ist=1:M
1134 for k=1:K
1135 if any(sum([QN(ist,k),UN(ist,k),TN(ist,k)])>0)
1136 JobClass{end+1,1} = self.sn{1}.classnames{k};
1137 Station{end+1,1} = self.sn{1}.nodenames{self.sn{1}.stationToNode(ist)};
1138 Qval(end+1) = QN(ist,k);
1139 Uval(end+1) = UN(ist,k);
1140 Tval(end+1) = TN(ist,k);
1141 end
1142 end
1143 end
1144 QLen = Qval(:); % we need to save first in a variable named like the column
1145 QT = Table(Station,JobClass,QLen);
1146 Util = Uval(:); % we need to save first in a variable named like the column
1147 UT = Table(Station,JobClass,Util);
1148 Tput = Tval(:); % we need to save first in a variable named like the column
1149 Station = categorical(Station);
1150 JobClass = categorical(JobClass);
1151 TT = Table(Station,JobClass,Tput);
1152 RespT = QLen ./ Tput;
1153 AvgTable = Table(Station,JobClass,QLen,Util,RespT,Tput);
1154 else
1155 Qval = zeros(M,K); Uval = zeros(M,K);
1156 JobClass = cell(K*M,1);
1157 Station = cell(K*M,1);
1158 for ist=1:M
1159 for k=1:K
1160 JobClass{(ist-1)*K+k} = Q{ist,k}.class.name;
1161 Station{(ist-1)*K+k} = Q{ist,k}.station.name;
1162 Qval((ist-1)*K+k) = QN(ist,k);
1163 Uval((ist-1)*K+k) = UN(ist,k);
1164 Tval((ist-1)*K+k) = TN(ist,k);
1165 end
1166 end
1167 Station = categorical(Station);
1168 JobClass = categorical(JobClass);
1169 QLen = Qval(:); % we need to save first in a variable named like the column
1170 QT = Table(Station,JobClass,QLen);
1171 Util = Uval(:); % we need to save first in a variable named like the column
1172 UT = Table(Station,JobClass,Util);
1173 Tput = Tval(:); % we need to save first in a variable named like the column
1174 TT = Table(Station,JobClass,Tput);
1175 RespT = QLen ./ Tput;
1176 AvgTable = Table(Station,JobClass,QLen,Util,RespT,Tput);
1177 end
1178 end
1179
1180 function envsn = getStruct(self)
1181 E = self.getNumberOfModels;
1182 envsn = cell(1,E);
1183 for e=1:E
1184 envsn{e} = self.ensemble{e}.getStruct;
1185 end
1186 end
1187
1188 function [T, segmentResults] = getSamplePathTable(self, samplePath)
1189 % [T, SEGMENTRESULTS] = GETSAMPLEPATHTABLE(SELF, SAMPLEPATH)
1190 %
1191 % Compute transient performance metrics for a sample path through
1192 % environment states. The method runs transient analysis for each
1193 % segment and extracts initial and final metric values.
1194 %
1195 % Inputs:
1196 % samplePath - Cell array where each row is {stage, duration}
1197 % stage: string (name) or integer (1-based index)
1198 % duration: positive scalar (time spent in stage)
1199 %
1200 % Outputs:
1201 % T - Table with columns: Segment, Stage, Duration, Station, JobClass,
1202 % InitQLen, InitUtil, InitTput, FinalQLen, FinalUtil, FinalTput
1203 % segmentResults - Cell array with detailed transient results per segment
1204 %
1205 % Example:
1206 % samplePath = {'Fast', 5.0; 'Slow', 10.0; 'Fast', 3.0};
1207 % [T, results] = solver.getSamplePathTable(samplePath);
1208
1209 % Validate input
1210 if isempty(samplePath)
1211 line_error(mfilename, 'Sample path cannot be empty.');
1212 end
1213
1214 % Initialize if needed
1215 if isempty(self.envObj.probEnv)
1216 self.init();
1217 end
1218
1219 E = self.getNumberOfModels;
1220 M = self.sn{1}.nstations;
1221 K = self.sn{1}.nclasses;
1222 nSegments = size(samplePath, 1);
1223 segmentResults = cell(nSegments, 1);
1224
1225 % Initialize queue lengths (uniform distribution for closed classes)
1226 Q_current = zeros(M, K);
1227 for k = 1:K
1228 if self.sn{1}.njobs(k) > 0 % closed class
1229 Q_current(:, k) = self.sn{1}.njobs(k) / M;
1230 end
1231 end
1232
1233 % Process each segment
1234 for seg = 1:nSegments
1235 % Resolve stage
1236 stageSpec = samplePath{seg, 1};
1237 duration = samplePath{seg, 2};
1238
1239 if ischar(stageSpec) || isstring(stageSpec)
1240 e = self.envObj.envGraph.findnode(stageSpec);
1241 if e == 0
1242 line_error(mfilename, sprintf('Stage "%s" not found.', stageSpec));
1243 end
1244 stageName = char(stageSpec);
1245 else
1246 e = stageSpec;
1247 if e < 1 || e > E
1248 line_error(mfilename, sprintf('Stage index %d out of range [1, %d].', e, E));
1249 end
1250 stageName = self.envObj.envGraph.Nodes.Name{e};
1251 end
1252
1253 if duration <= 0
1254 line_error(mfilename, 'Duration must be positive.');
1255 end
1256
1257 % Initialize model from current queue lengths
1258 self.ensemble{e}.initFromMarginal(Q_current);
1259
1260 % Set solver timespan
1261 self.solvers{e}.options.timespan = [0, duration];
1262 self.solvers{e}.reset();
1263
1264 % Run transient analysis
1265 [Qt, Ut, Tt] = self.ensemble{e}.getTranHandles;
1266 [QNt, UNt, TNt] = self.solvers{e}.getTranAvg(Qt, Ut, Tt);
1267
1268 % Extract initial and final metrics
1269 initQ = zeros(M, K);
1270 initU = zeros(M, K);
1271 initT = zeros(M, K);
1272 finalQ = zeros(M, K);
1273 finalU = zeros(M, K);
1274 finalT = zeros(M, K);
1275
1276 for i = 1:M
1277 for r = 1:K
1278 Qir = QNt{i, r};
1279 Uir = UNt{i, r};
1280 Tir = TNt{i, r};
1281
1282 if isstruct(Qir) && isfield(Qir, 'metric') && ~isempty(Qir.metric)
1283 initQ(i, r) = Qir.metric(1);
1284 finalQ(i, r) = Qir.metric(end);
1285 end
1286 if isstruct(Uir) && isfield(Uir, 'metric') && ~isempty(Uir.metric)
1287 initU(i, r) = Uir.metric(1);
1288 finalU(i, r) = Uir.metric(end);
1289 end
1290 if isstruct(Tir) && isfield(Tir, 'metric') && ~isempty(Tir.metric)
1291 initT(i, r) = Tir.metric(1);
1292 finalT(i, r) = Tir.metric(end);
1293 end
1294 end
1295 end
1296
1297 % Store segment results
1298 segmentResults{seg}.stage = e;
1299 segmentResults{seg}.stageName = stageName;
1300 segmentResults{seg}.duration = duration;
1301 segmentResults{seg}.QNt = QNt;
1302 segmentResults{seg}.UNt = UNt;
1303 segmentResults{seg}.TNt = TNt;
1304 segmentResults{seg}.initQ = initQ;
1305 segmentResults{seg}.initU = initU;
1306 segmentResults{seg}.initT = initT;
1307 segmentResults{seg}.finalQ = finalQ;
1308 segmentResults{seg}.finalU = finalU;
1309 segmentResults{seg}.finalT = finalT;
1310
1311 % Update current queue lengths for next segment
1312 Q_current = finalQ;
1313 end
1314
1315 % Build output table
1316 Segment = [];
1317 Stage = {};
1318 Duration = [];
1319 Station = {};
1320 JobClass = {};
1321 InitQLen = [];
1322 InitUtil = [];
1323 InitTput = [];
1324 FinalQLen = [];
1325 FinalUtil = [];
1326 FinalTput = [];
1327
1328 for seg = 1:nSegments
1329 res = segmentResults{seg};
1330 for i = 1:M
1331 for r = 1:K
1332 % Only include rows with non-zero metrics
1333 if any([res.initQ(i,r), res.initU(i,r), res.initT(i,r), ...
1334 res.finalQ(i,r), res.finalU(i,r), res.finalT(i,r)] > 0)
1335 Segment(end+1, 1) = seg;
1336 Stage{end+1, 1} = res.stageName;
1337 Duration(end+1, 1) = res.duration;
1338 Station{end+1, 1} = self.sn{1}.nodenames{self.sn{1}.stationToNode(i)};
1339 JobClass{end+1, 1} = self.sn{1}.classnames{r};
1340 InitQLen(end+1, 1) = res.initQ(i, r);
1341 InitUtil(end+1, 1) = res.initU(i, r);
1342 InitTput(end+1, 1) = res.initT(i, r);
1343 FinalQLen(end+1, 1) = res.finalQ(i, r);
1344 FinalUtil(end+1, 1) = res.finalU(i, r);
1345 FinalTput(end+1, 1) = res.finalT(i, r);
1346 end
1347 end
1348 end
1349 end
1350
1351 Stage = categorical(Stage);
1352 Station = categorical(Station);
1353 JobClass = categorical(JobClass);
1354 T = Table(Segment, Stage, Duration, Station, JobClass, ...
1355 InitQLen, InitUtil, InitTput, FinalQLen, FinalUtil, FinalTput);
1356 end
1357 function [allMethods] = listValidMethods(self)
1358 % allMethods = LISTVALIDMETHODS()
1359 % List valid methods for this solver
1360 sn = self.model.getStruct();
1361 allMethods = {'default'};
1362 end
1363
1364 function Q = roundMarginalForDiscreteSolver(self, Q, snRef)
1365 % ROUNDMARGINALFORDISCRETESOLVER Round fractional queue lengths
1366 % to integers using the largest remainder method, preserving
1367 % closed chain populations exactly.
1368 for c = 1:snRef.nchains
1369 chainClasses = find(snRef.chains(c,:) > 0);
1370 njobs_chain = sum(snRef.njobs(chainClasses));
1371 if isinf(njobs_chain)
1372 % Open chain: simple rounding
1373 for idx = 1:length(chainClasses)
1374 k = chainClasses(idx);
1375 for i = 1:size(Q,1)
1376 Q(i,k) = round(Q(i,k));
1377 end
1378 end
1379 else
1380 % Closed chain: largest remainder method
1381 vals = [];
1382 indices = [];
1383 for i = 1:size(Q,1)
1384 for idx = 1:length(chainClasses)
1385 k = chainClasses(idx);
1386 vals(end+1) = Q(i,k);
1387 indices(end+1,:) = [i, k];
1388 end
1389 end
1390 floored = floor(vals);
1391 remainders = vals - floored;
1392 deficit = njobs_chain - sum(floored);
1393 deficit = round(deficit); % ensure integer
1394 if deficit > 0
1395 [~, sortIdx] = sort(remainders, 'descend');
1396 for d = 1:min(deficit, length(sortIdx))
1397 floored(sortIdx(d)) = floored(sortIdx(d)) + 1;
1398 end
1399 end
1400 for j = 1:length(vals)
1401 Q(indices(j,1), indices(j,2)) = floored(j);
1402 end
1403 end
1404 end
1405 end
1406 end
1407
1408 methods (Static)
1409
1410 function [bool, featSupported] = supports(model)
1411 % [BOOL, FEATSUPPORTED] = SUPPORTS(MODEL)
1412
1413 featUsed = model.getUsedLangFeatures();
1414
1415 featSupported = SolverFeatureSet;
1416
1417 % Nodes
1418 featSupported.setTrue('ClassSwitch');
1419 featSupported.setTrue('Delay');
1420 featSupported.setTrue('DelayStation');
1421 featSupported.setTrue('Queue');
1422 featSupported.setTrue('Sink');
1423 featSupported.setTrue('Source');
1424
1425 % Distributions
1426 featSupported.setTrue('Coxian');
1427 featSupported.setTrue('Cox2');
1428 featSupported.setTrue('Erlang');
1429 featSupported.setTrue('Exp');
1430 featSupported.setTrue('HyperExp');
1431
1432 % Sections
1433 featSupported.setTrue('StatelessClassSwitcher'); % Section
1434 featSupported.setTrue('InfiniteServer'); % Section
1435 featSupported.setTrue('SharedServer'); % Section
1436 featSupported.setTrue('Buffer'); % Section
1437 featSupported.setTrue('Dispatcher'); % Section
1438 featSupported.setTrue('Server'); % Section (Non-preemptive)
1439 featSupported.setTrue('JobSink'); % Section
1440 featSupported.setTrue('RandomSource'); % Section
1441 featSupported.setTrue('ServiceTunnel'); % Section
1442
1443 % Scheduling strategy
1444 featSupported.setTrue('SchedStrategy_INF');
1445 featSupported.setTrue('SchedStrategy_PS');
1446 featSupported.setTrue('SchedStrategy_FCFS');
1447 featSupported.setTrue('RoutingStrategy_PROB');
1448 featSupported.setTrue('RoutingStrategy_RAND');
1449 featSupported.setTrue('RoutingStrategy_RROBIN'); % with SolverJMT
1450
1451 % Customer Classes
1452 featSupported.setTrue('ClosedClass');
1453 featSupported.setTrue('OpenClass');
1454
1455 bool = SolverFeatureSet.supports(featSupported, featUsed);
1456 end
1457 end
1458
1459 methods (Static)
1460 % ensemble solver options
1461 function options = defaultOptions()
1462 % OPTIONS = DEFAULTOPTIONS()
1463 options = SolverOptions('ENV');
1464 end
1465
1466 function libs = getLibrariesUsed(sn, options)
1467 % GETLIBRARIESUSED Get list of external libraries used by ENV solver
1468 % ENV uses internal algorithms, no external library attribution needed
1469 libs = {};
1470 end
1471 end
1472end
1473