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