LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
SolverLN.m
1classdef SolverLN < EnsembleSolver
2 % SolverLN Layered network solver for hierarchical performance models
3 %
4 % SolverLN implements analysis of layered queueing networks (LQNs) which model
5 % hierarchical software systems with clients, application servers, and resource
6 % layers. It uses iterative decomposition to solve multi-layer models by
7 % analyzing each layer separately and propagating service demands between layers.
8 %
9 % @brief Layered network solver for hierarchical software performance models
10 %
11 % Example:
12 % @code
13 % solver = SolverLN(layered_model, 'maxIter', 100);
14 % solver.runAnalyzer(); % Iterative layer analysis
15 % metrics = solver.getEnsembleAvg(); % Layer performance metrics
16 % @endcode
17 %
18 % Copyright (c) 2012-2026, Imperial College London
19 % All rights reserved.
20
21 properties %(Hidden) % registries of quantities to update at every iteration
22 nlayers; % number of model layers
23 lqn; % lqn data structure
24 hasconverged; % true if last iteration converged, false otherwise
25 averagingstart; % iteration at which result averaging started
26 idxhash; % ensemble model associated to host or task
27 servtmatrix; % auxiliary matrix to determine entry servt
28 ptaskcallers; % probability that a task is called by a given task, directly or indirectly (remotely)
29 ptaskcallers_step; % probability that a task is called by a given task, directly or indirectly (remotely) up to a given step distance
30 ilscaling; % interlock scalings
31 njobs; % number of jobs for each caller in a given submodel
32 njobsorig; % number of jobs for each caller at layer build time
33 routereset; % models that require hard reset of service chains
34 svcreset; % models that require hard reset of service process
35 maxitererr; % maximum error at current iteration over all layers
36 % Under-relaxation state for convergence improvement
37 relax_omega; % Current relaxation factor
38 relax_err_history; % Error history for adaptive mode
39 servt_prev; % Previous service times for relaxation
40 residt_prev; % Previous residence times for relaxation
41 tput_prev; % Previous throughputs for relaxation
42 thinkt_prev; % Previous think times for relaxation
43 callservt_prev; % Previous call service times for relaxation
44 % MOL (Method of Layers) properties for hierarchical iteration
45 hostLayerIndices; % Indices of host (processor) layers in ensemble
46 taskLayerIndices; % Indices of task layers in ensemble
47 util_prev_host; % Previous processor utilizations (for delta computation)
48 util_prev_task; % Previous task utilizations (for delta computation)
49 mol_it_host_outer; % MOL host outer iteration counter
50 mol_it_task_inner; % MOL task inner iteration counter
51 % Phase-2 support properties
52 hasPhase2; % Flag: model has phase-2 activities
53 servt_ph1; % Phase-1 service time per activity (nidx x 1)
54 servt_ph2; % Phase-2 service time per activity (nidx x 1)
55 util_ph1; % Phase-1 utilization per entry
56 util_ph2; % Phase-2 utilization per entry
57 prOvertake; % Overtaking probability per entry (nentries x 1)
58 end
59
60 properties %(Hidden) % performance metrics and related processes
61 util;
62 util_ilock; % interlock matrix (ntask x ntask), element (i,j) says how much the utilization of task i is imputed to task j
63 tput;
64 tputproc;
65 servt; % this is the mean service time of an activity, which is the response time at the lower layer (if applicable)
66 residt; % this is the residence time at the lower layer (if applicable)
67 servtproc; % this is the service time process with mean fitted to the servt value
68 servtcdf; % this is the cdf of the service time process
69 thinkt;
70 thinkproc;
71 thinktproc;
72 entryproc;
73 entrycdfrespt;
74 callresidt;
75 callservt;
76 callservtproc;
77 callservtcdf;
78 ignore; % elements to be ignored (e.g., components disconnected from a REF node)
79 end
80
81 properties %(Access = protected, Hidden) % registries of quantities to update at every iteration
82 arvproc_classes_updmap; % [modelidx, actidx, node, class]
83 thinkt_classes_updmap; % [modelidx, actidx, node, class]
84 actthinkt_classes_updmap; % [modelidx, actidx, node, class] for activity think-times
85 servt_classes_updmap; % [modelidx, actidx, node, class]
86 call_classes_updmap; % [modelidx, callidx, node, class]
87 route_prob_updmap; % [modelidx, actidxfrom, actidxto, nodefrom, nodeto, classfrom, classto]
88 unique_route_prob_updmap; % auxiliary cache of unique route_prob_updmap rows
89 solverFactory; % function handle to create layer solvers
90 end
91
92 methods
93 function self = SolverLN(lqnmodel, solverFactory, varargin)
94 % SELF = SOLVERLN(MODEL,SOLVERFACTORY,VARARGIN)
95 self@EnsembleSolver(lqnmodel, mfilename);
96
97 if any(cellfun(@(s) strcmpi(s,'java'), varargin))
98 self.obj = JLINE.SolverLN(JLINE.from_line_layered_network(lqnmodel));
99 self.obj.options.verbose = jline.VerboseLevel.SILENT;
100 else
101 % Default solver factory: Use JMT for open networks, MVA for closed networks
102 defaultSolverFactory = @(m) adaptiveSolverFactory(m, self.options);
103
104 if nargin == 1 %case SolverLN(model)
105 solverFactory = defaultSolverFactory;
106 self.setOptions(SolverLN.defaultOptions);
107 elseif nargin>1 && isstruct(solverFactory)
108 options = solverFactory;
109 self.setOptions(options);
110 solverFactory = defaultSolverFactory;
111 elseif nargin>2 % case SolverLN(model,'opt1',...)
112 if ischar(solverFactory)
113 inputvar = {solverFactory,varargin{:}}; %#ok<CCAT>
114 solverFactory = defaultSolverFactory;
115 else % case SolverLN(model, solverFactory, 'opt1',...)
116 inputvar = varargin;
117 end
118 self.setOptions(Solver.parseOptions(inputvar, SolverLN.defaultOptions));
119 else %case SolverLN(model,solverFactory)
120 self.setOptions(SolverLN.defaultOptions);
121 end
122 self.lqn = lqnmodel.getStruct();
123
124 % Detect and initialize phase-2 support
125 if isfield(self.lqn, 'actphase') && any(self.lqn.actphase > 1)
126 self.hasPhase2 = true;
127 self.servt_ph1 = zeros(self.lqn.nidx, 1);
128 self.servt_ph2 = zeros(self.lqn.nidx, 1);
129 self.util_ph1 = zeros(self.lqn.nidx, 1);
130 self.util_ph2 = zeros(self.lqn.nidx, 1);
131 self.prOvertake = zeros(self.lqn.nentries, 1);
132 else
133 self.hasPhase2 = false;
134 end
135
136 self.construct();
137 for e=1:self.getNumberOfModels
138 if numel(find(self.lqn.isfunction == 1))
139 if ~isempty(self.ensemble{e}.stations{2}.setupTime)
140 solverFactory = @(m) SolverMAM(m,'verbose',false,'method','dec.poisson');
141 else
142 solverFactory = @(m) SolverAuto(m,'verbose',false);
143 end
144 self.setSolver(solverFactory(self.ensemble{e}),e);
145 else
146 self.setSolver(solverFactory(self.ensemble{e}),e);
147 end
148 end
149 self.solverFactory = solverFactory; % Store for later use
150 end
151 end
152
153 function runtime = runAnalyzer(self, options) %#ok<INUSD> % generic method to run the solver
154 line_error(mfilename,'Use getEnsembleAvg instead.');
155 end
156
157 function sn = getStruct(self)
158 % SN = GETSTRUCT()
159
160 % Get data structure summarizing the model
161 sn = self.model.getStruct();
162 end
163
164 function construct(self)
165 % mark down to ignore unreachable disconnected components
166 self.ignore = false(self.lqn.nidx,1);
167 [~,wccs] = weaklyconncomp(self.lqn.graph'+self.lqn.graph);
168 uwccs = unique(wccs);
169 if length(uwccs)>1
170 % the model has disconnected submodels
171 wccref = false(1,length(uwccs));
172 for t=1:self.lqn.ntasks
173 tidx = self.lqn.tshift+t;
174 if self.lqn.sched(tidx) == SchedStrategy.REF
175 wccref(wccs(tidx)) = true;
176 end
177 end
178 if any(wccref==false)
179 for dw=find(wccref==false) % disconnected component
180 self.ignore(find(wccs==dw)) = true;
181 end
182 end
183 end
184
185 % initialize internal data structures
186 self.entrycdfrespt = cell(length(self.lqn.nentries),1);
187 self.hasconverged = false;
188
189 % initialize svc and think times
190 self.servtproc = self.lqn.hostdem;
191 self.thinkproc = self.lqn.think;
192 self.callservtproc = cell(self.lqn.ncalls,1);
193 for cidx = 1:self.lqn.ncalls
194 self.callservtproc{cidx} = self.lqn.hostdem{self.lqn.callpair(cidx,2)};
195 end
196
197 % perform layering
198 self.njobs = zeros(self.lqn.tshift + self.lqn.ntasks, self.lqn.tshift + self.lqn.ntasks);
199 buildLayers(self); % build layers
200 self.njobsorig = self.njobs;
201 self.nlayers = length(self.ensemble);
202
203 % initialize data structures for interlock correction
204 self.ptaskcallers = zeros(self.lqn.nhosts+self.lqn.ntasks, self.lqn.nhosts+self.lqn.ntasks);
205 self.ptaskcallers_step = cell(1,self.nlayers+1);
206 for step=1:self.nlayers % upper bound on maximum dag height
207 self.ptaskcallers_step{step} = zeros(self.lqn.nhosts+self.lqn.ntasks, self.lqn.nhosts+self.lqn.ntasks);
208 end
209
210 % layering generates update maps that we use here to cache the elements that need reset
211 self.routereset = unique(self.idxhash(self.route_prob_updmap(:,1)))';
212 self.svcreset = unique(self.idxhash(self.thinkt_classes_updmap(:,1)))';
213 self.svcreset = union(self.svcreset,unique(self.idxhash(self.call_classes_updmap(:,1)))');
214 end
215
216 function self = reset(self)
217 % no-op
218 end
219
220 bool = converged(self, it); % convergence test at iteration it
221
222 function init(self) % operations before starting to iterate
223 % INIT() % OPERATIONS BEFORE STARTING TO ITERATE
224 self.unique_route_prob_updmap = unique(self.route_prob_updmap(:,1))';
225 self.tput = zeros(self.lqn.nidx,1);
226 self.tputproc = cell(self.lqn.nidx,1);
227 self.util = zeros(self.lqn.nidx,1);
228 self.servt = zeros(self.lqn.nidx,1);
229 self.servtmatrix = getEntryServiceMatrix(self);
230
231 for e= 1:self.nlayers
232 self.solvers{e}.enableChecks=false;
233 end
234
235 % Initialize under-relaxation state
236 relax_mode = self.options.config.relax;
237 switch relax_mode
238 case {'auto'}
239 self.relax_omega = 1.0; % Start without relaxation
240 case {'fixed', 'adaptive'}
241 self.relax_omega = self.options.config.relax_factor;
242 otherwise % 'none' or unrecognized
243 self.relax_omega = 1.0; % No relaxation
244 end
245 self.relax_err_history = [];
246 self.servt_prev = NaN(self.lqn.nidx, 1);
247 self.residt_prev = NaN(self.lqn.nidx, 1);
248 self.tput_prev = NaN(self.lqn.nidx, 1);
249 self.thinkt_prev = NaN(self.lqn.nidx, 1);
250 self.callservt_prev = NaN(self.lqn.ncalls, 1);
251
252 % Initialize MOL-specific state
253 self.mol_it_host_outer = 0;
254 self.mol_it_task_inner = 0;
255 self.util_prev_host = zeros(self.lqn.nhosts, 1);
256 self.util_prev_task = zeros(self.lqn.ntasks, 1);
257 end
258
259
260 function pre(self, it) %#ok<INUSD> % operations before an iteration
261 % PRE(IT) % OPERATIONS BEFORE AN ITERATION
262 % no-op
263 end
264
265 function [result, runtime] = analyze(self, it, e) %#ok<INUSD>
266 % [RESULT, RUNTIME] = ANALYZE(IT, E)
267 T0 = tic;
268 result = struct();
269 %jresult = struct();
270 if e==1 && self.solvers{e}.options.verbose
271 line_printf('\n');
272 end
273
274 % Protection for unstable queues during LN iterations
275 % If a solver fails (e.g., due to queue instability with open arrivals),
276 % use results from previous iteration if available and continue
277 try
278 [result.QN, result.UN, result.RN, result.TN, result.AN, result.WN] = self.solvers{e}.getAvg();
279 catch ME
280 if it > 1 && ~isempty(self.results) && size(self.results, 1) >= (it-1) && size(self.results, 2) >= e
281 if self.solvers{e}.options.verbose
282 warning('LINE:SolverLN:Instability', ...
283 'Layer %d at iteration %d encountered instability (possibly due to high service demand with open arrivals). Using previous iteration values and continuing.', ...
284 e, it);
285 end
286 % Use results from previous iteration
287 prevResult = self.results{it-1, e};
288 result.QN = prevResult.QN;
289 result.UN = prevResult.UN;
290 result.RN = prevResult.RN;
291 result.TN = prevResult.TN;
292 result.AN = prevResult.AN;
293 result.WN = prevResult.WN;
294 else
295 % First iteration or no previous results, re-throw the exception
296 error('LINE:SolverLN:FirstIterationFailure', ...
297 'Layer %d failed at iteration %d with no previous iteration to fall back on: %s', ...
298 e, it, ME.message);
299 end
300 end
301 runtime = toc(T0);
302 end
303
304 function post(self, it) % operations after an iteration
305 % POST(IT) % OPERATIONS AFTER AN ITERATION
306 % convert the results of QNs into layer metrics
307
308 self.updateMetrics(it);
309
310 % recompute think times
311 self.updateThinkTimes(it);
312
313 if self.options.config.interlocking
314 % recompute layer populations
315 self.updatePopulations(it);
316 end
317
318 % update the model parameters
319 self.updateLayers(it);
320
321 % update entry selection and cache routing probabilities within callers
322 self.updateRoutingProbabilities(it);
323
324 % reset all layers with routing probability changes
325 for e= self.routereset
326 self.ensemble{e}.refreshChains();
327 self.solvers{e}.reset();
328 end
329
330 % refresh visits and network model parameters
331 for e= self.svcreset
332 switch self.solvers{e}.name
333 case {'SolverMVA', 'SolverNC'} %leaner than refreshProcesses, no need to refresh phases
334 % note: this does not refresh the sn.proc field, only sn.rates and sn.scv
335 switch self.options.method
336 case 'default'
337 self.ensemble{e}.refreshRates();
338 case 'moment3'
339 self.ensemble{e}.refreshProcesses();
340 end
341 otherwise
342 self.ensemble{e}.refreshProcesses();
343 end
344 self.solvers{e}.reset(); % commenting this out des not seem to produce a problem, but it goes faster with it
345 end
346
347 % this is required to handle population changes due to interlocking
348 if self.options.config.interlocking
349 for e=1:self.nlayers
350 self.ensemble{e}.refreshJobs();
351 end
352 end
353
354 if it==1
355 % now disable all solver support checks for future iterations
356 for e=1:length(self.ensemble)
357 self.solvers{e}.setChecks(false);
358 end
359 end
360 end
361
362
363 function finish(self) % operations after iterations are completed
364 % FINISH() % OPERATIONS AFTER INTERATIONS ARE COMPLETED
365 E = size(self.results,2);
366 for e=1:E
367 s = self.solvers{e};
368 s.getAvg();
369 self.solvers{e} = s;
370 end
371 self.model.ensemble = self.ensemble;
372 end
373
374 function [QNlqn_t, UNlqn_t, TNlqn_t] = getTranAvg(self)
375 self.getAvg;
376 QNclass_t = {};
377 UNclass_t = {};
378 TNclass_t = {};
379 QNlqn_t = cell(0,0);
380 for e=1:self.nlayers
381 [crows, ccols] = size(QNlqn_t);
382 [QNclass_t{e}, UNclass_t{e}, TNclass_t{e}] = self.solvers{e}.getTranAvg();
383 QNlqn_t(crows+1:crows+size(QNclass_t{e},1),ccols+1:ccols+size(QNclass_t{e},2)) = QNclass_t{e};
384 UNlqn_t(crows+1:crows+size(UNclass_t{e},1),ccols+1:ccols+size(UNclass_t{e},2)) = UNclass_t{e};
385 TNlqn_t(crows+1:crows+size(TNclass_t{e},1),ccols+1:ccols+size(TNclass_t{e},2)) = TNclass_t{e};
386 end
387 end
388
389 function varargout = getAvg(varargin)
390 % [QN,UN,RN,TN,AN,WN] = GETAVG(SELF,~,~,~,~,USELQNSNAMING)
391 [varargout{1:nargout}] = getEnsembleAvg( varargin{:} );
392 end
393
394 function [cdfRespT] = getCdfRespT(self)
395 if isempty(self.entrycdfrespt{1})
396 % save user-specified method to temporary variable
397 curMethod = self.getOptions.method;
398 % run with moment 3
399 self.options.method = 'moment3';
400 self.getAvg();
401 % restore user-specified method
402 self.options.method = curMethod;
403 end
404 cdfRespT = self.entrycdfrespt;
405 end
406
407 function [AvgTable,QT,UT,RT,WT,AT,TT] = getAvgTable(self)
408 % [AVGTABLE,QT,UT,RT,WT,TT] = GETAVGTABLE(USELQNSNAMING)
409 if (GlobalConstants.DummyMode)
410 [AvgTable, QT, UT, RT, TT, WT] = deal([]);
411 return
412 end
413
414 if ~isempty(self.obj)
415 avgTable = self.obj.getEnsembleAvg();
416 [QN,UN,RN,WN,AN,TN] = JLINE.arrayListToResults(avgTable);
417 else
418 [QN,UN,RN,TN,AN,WN] = getAvg(self);
419 end
420
421 % attempt to sanitize small numerical perturbations
422 variables = {QN, UN, RN, TN, AN, WN}; % Put all variables in a cell array
423 for i = 1:length(variables)
424 rVar = round(variables{i} * 10);
425 toRound = abs(variables{i} * 10 - rVar) < GlobalConstants.CoarseTol * variables{i} * 10;
426 variables{i}(toRound) = rVar(toRound) / 10;
427 variables{i}(variables{i}<=GlobalConstants.FineTol) = 0;
428 end
429 [QN, UN, RN, TN, AN, WN] = deal(variables{:}); % Assign the modified values back to the original variables
430
431 %%
432 Node = label(self.lqn.names);
433 O = length(Node);
434 NodeType = label(O,1);
435 for o = 1:O
436 switch self.lqn.type(o)
437 case LayeredNetworkElement.PROCESSOR
438 NodeType(o,1) = label({'Processor'});
439 case LayeredNetworkElement.TASK
440 if self.lqn.isref(o)
441 NodeType(o,1) = label({'RefTask'});
442 else
443 NodeType(o,1) = label({'Task'});
444 end
445 case LayeredNetworkElement.ENTRY
446 NodeType(o,1) = label({'Entry'});
447 case LayeredNetworkElement.ACTIVITY
448 NodeType(o,1) = label({'Activity'});
449 case LayeredNetworkElement.CALL
450 NodeType(o,1) = label({'Call'});
451 end
452 end
453 QLen = QN;
454 QT = Table(Node,QLen);
455 Util = UN;
456 UT = Table(Node,Util);
457 RespT = RN;
458 RT = Table(Node,RespT);
459 Tput = TN;
460 TT = Table(Node,Tput);
461 %SvcT = SN;
462 %ST = Table(Node,SvcT);
463 %ProcUtil = PN;
464 %PT = Table(Node,ProcUtil);
465 ResidT = WN;
466 WT = Table(Node,ResidT);
467 ArvR = AN;
468 AT = Table(Node,ArvR);
469 AvgTable = Table(Node, NodeType, QLen, Util, RespT, ResidT, ArvR, Tput);%, ProcUtil, SvcT);
470 end
471
472 function [AvgTable,QT,UT,RT,WT,AT,TT] = avgTable(self)
473 % AVGTABLE Alias for getAvgTable
474 [AvgTable,QT,UT,RT,WT,AT,TT] = self.getAvgTable();
475 end
476
477 function [AvgTable,QT,UT,RT,WT,AT,TT] = avgT(self)
478 % AVGT Short alias for getAvgTable
479 [AvgTable,QT,UT,RT,WT,AT,TT] = self.getAvgTable();
480 end
481
482 function [AvgTable,QT,UT,RT,WT,AT,TT] = aT(self)
483 % AT Short alias for getAvgTable (MATLAB-compatible)
484 [AvgTable,QT,UT,RT,WT,AT,TT] = self.getAvgTable();
485 end
486 end
487
488 methods
489 [QN,UN,RN,TN,AN,WN] = getEnsembleAvg(self);
490
491 function [bool, featSupported] = supports(model)
492 % [BOOL, FEATSUPPORTED] = SUPPORTS(MODEL)
493 % This method cannot be static as otherwise it cannot access self.solvers{e}
494 ensemble = model.getEnsemble;
495 featSupported = cell(length(ensemble),1);
496 bool = true;
497 for e = 1:length(ensemble)
498 [solverSupports,featSupported{e}] = self.solvers{e}.supports(ensemble{e});
499 bool = bool && solverSupports;
500 end
501 end
502 end
503
504 methods (Hidden)
505 buildLayers(self, lqn, resptproc, callservtproc);
506 buildLayersRecursive(self, idx, callers, ishostlayer);
507 updateLayers(self, it);
508 updatePopulations(self, it);
509 updateThinkTimes(self, it);
510 updateMetrics(self, it);
511 updateRoutingProbabilities(self, it);
512 svcmatrix = getEntryServiceMatrix(self)
513 prOt = overtake_prob(self, eidx); % Phase-2 overtaking probability
514
515 function delta = computeTaskDelta(self)
516 % COMPUTETASKDELTA Compute max queue-length change for task layers (MOL inner loop)
517 %
518 % Returns the maximum queue-length change across task layers,
519 % normalized by total jobs, similar to converged.m logic.
520 results = self.results;
521 it = size(results, 1);
522 if it < 2
523 delta = Inf;
524 return;
525 end
526
527 delta = 0;
528 for e = self.taskLayerIndices
529 metric = results{it, e}.QN;
530 metric_1 = results{it-1, e}.QN;
531 N = sum(self.ensemble{e}.getNumberOfJobs);
532 if N > 0
533 try
534 iterErr = max(abs(metric(:) - metric_1(:))) / N;
535 catch
536 iterErr = 0;
537 end
538 delta = max(delta, iterErr);
539 end
540 end
541 end
542
543 function delta = computeHostDelta(self)
544 % COMPUTEHOSTDELTA Compute max queue-length change for host layers (MOL outer loop)
545 %
546 % Returns the maximum queue-length change across host layers,
547 % normalized by total jobs, similar to converged.m logic.
548 results = self.results;
549 it = size(results, 1);
550 if it < 2
551 delta = Inf;
552 return;
553 end
554
555 delta = 0;
556 for e = self.hostLayerIndices
557 metric = results{it, e}.QN;
558 metric_1 = results{it-1, e}.QN;
559 N = sum(self.ensemble{e}.getNumberOfJobs);
560 if N > 0
561 try
562 iterErr = max(abs(metric(:) - metric_1(:))) / N;
563 catch
564 iterErr = 0;
565 end
566 delta = max(delta, iterErr);
567 end
568 end
569 end
570
571 function postTaskLayers(self, it)
572 % POSTTASKLAYERS Post-iteration updates for task layers only (MOL inner loop)
573 %
574 % Updates think times, layer parameters, and routing probabilities
575 % after task layer analysis, then resets task layer solvers.
576
577 self.updateThinkTimes(it);
578 if self.options.config.interlocking
579 self.updatePopulations(it);
580 end
581 self.updateLayers(it);
582 self.updateRoutingProbabilities(it);
583
584 % Reset task layer solvers
585 for e = self.taskLayerIndices
586 switch self.solvers{e}.name
587 case {'SolverMVA', 'SolverNC'}
588 switch self.options.method
589 case 'mol'
590 self.ensemble{e}.refreshRates();
591 otherwise
592 self.ensemble{e}.refreshRates();
593 end
594 otherwise
595 self.ensemble{e}.refreshProcesses();
596 end
597 self.solvers{e}.reset();
598 end
599
600 if self.options.config.interlocking
601 for e = self.taskLayerIndices
602 self.ensemble{e}.refreshJobs();
603 end
604 end
605 end
606
607 function postHostLayers(self, it)
608 % POSTHOSTLAYERS Post-iteration updates for host layers (MOL outer loop)
609 %
610 % Updates think times, layer parameters, and routing probabilities
611 % after host layer analysis, then resets host layer solvers.
612
613 self.updateThinkTimes(it);
614 if self.options.config.interlocking
615 self.updatePopulations(it);
616 end
617 self.updateLayers(it);
618 self.updateRoutingProbabilities(it);
619
620 % Reset host layer solvers
621 for e = self.hostLayerIndices
622 switch self.solvers{e}.name
623 case {'SolverMVA', 'SolverNC'}
624 self.ensemble{e}.refreshRates();
625 otherwise
626 self.ensemble{e}.refreshProcesses();
627 end
628 self.solvers{e}.reset();
629 end
630
631 if self.options.config.interlocking
632 for e = self.hostLayerIndices
633 self.ensemble{e}.refreshJobs();
634 end
635 end
636 end
637 end
638
639 methods
640 function state = get_state(self)
641 % GET_STATE Export current solver state for continuation
642 %
643 % STATE = GET_STATE() returns a struct containing the current
644 % solution state, which can be used to continue iteration with
645 % a different solver via set_state().
646 %
647 % The exported state includes:
648 % - Service time processes (servtproc)
649 % - Think time processes (thinktproc)
650 % - Call service time processes (callservtproc)
651 % - Throughput processes (tputproc)
652 % - Performance metrics (util, tput, servt, residt, etc.)
653 % - Relaxation state
654 % - Last iteration results
655 %
656 % Example:
657 % solver1 = SolverLN(model, @(m) SolverMVA(m));
658 % solver1.getEnsembleAvg();
659 % state = solver1.get_state();
660 %
661 % solver2 = SolverLN(model, @(m) SolverJMT(m));
662 % solver2.set_state(state);
663 % solver2.getEnsembleAvg(); % Continues from MVA solution
664
665 state = struct();
666
667 % Service/think time processes
668 state.servtproc = self.servtproc;
669 state.thinktproc = self.thinktproc;
670 state.callservtproc = self.callservtproc;
671 state.tputproc = self.tputproc;
672 state.entryproc = self.entryproc;
673
674 % Performance metrics
675 state.util = self.util;
676 state.tput = self.tput;
677 state.servt = self.servt;
678 state.residt = self.residt;
679 state.thinkt = self.thinkt;
680 state.callresidt = self.callresidt;
681 state.callservt = self.callservt;
682
683 % Relaxation state
684 state.relax_omega = self.relax_omega;
685 state.servt_prev = self.servt_prev;
686 state.residt_prev = self.residt_prev;
687 state.tput_prev = self.tput_prev;
688 state.thinkt_prev = self.thinkt_prev;
689 state.callservt_prev = self.callservt_prev;
690
691 % Results from last iteration
692 state.results = self.results;
693
694 % Interlock data
695 state.njobs = self.njobs;
696 state.ptaskcallers = self.ptaskcallers;
697 state.ilscaling = self.ilscaling;
698 end
699
700 function set_state(self, state)
701 % SET_STATE Import solution state for continuation
702 %
703 % SET_STATE(STATE) initializes the solver with a previously
704 % exported state, allowing iteration to continue from where
705 % a previous solver left off.
706 %
707 % This enables hybrid solving schemes where fast solvers (MVA)
708 % provide initial estimates and accurate solvers (JMT, DES)
709 % refine the solution.
710 %
711 % Example:
712 % solver1 = SolverLN(model, @(m) SolverMVA(m));
713 % solver1.getEnsembleAvg();
714 % state = solver1.get_state();
715 %
716 % solver2 = SolverLN(model, @(m) SolverJMT(m));
717 % solver2.set_state(state);
718 % solver2.getEnsembleAvg(); % Continues from MVA solution
719
720 % Service/think time processes
721 self.servtproc = state.servtproc;
722 self.thinktproc = state.thinktproc;
723 self.callservtproc = state.callservtproc;
724 self.tputproc = state.tputproc;
725 if isfield(state, 'entryproc')
726 self.entryproc = state.entryproc;
727 end
728
729 % Performance metrics
730 self.util = state.util;
731 self.tput = state.tput;
732 self.servt = state.servt;
733 if isfield(state, 'residt')
734 self.residt = state.residt;
735 end
736 if isfield(state, 'thinkt')
737 self.thinkt = state.thinkt;
738 end
739 if isfield(state, 'callresidt')
740 self.callresidt = state.callresidt;
741 end
742 if isfield(state, 'callservt')
743 self.callservt = state.callservt;
744 end
745
746 % Relaxation state
747 if isfield(state, 'relax_omega')
748 self.relax_omega = state.relax_omega;
749 end
750 if isfield(state, 'servt_prev')
751 self.servt_prev = state.servt_prev;
752 end
753 if isfield(state, 'residt_prev')
754 self.residt_prev = state.residt_prev;
755 end
756 if isfield(state, 'tput_prev')
757 self.tput_prev = state.tput_prev;
758 end
759 if isfield(state, 'thinkt_prev')
760 self.thinkt_prev = state.thinkt_prev;
761 end
762 if isfield(state, 'callservt_prev')
763 self.callservt_prev = state.callservt_prev;
764 end
765
766 % Results
767 if isfield(state, 'results')
768 self.results = state.results;
769 end
770
771 % Interlock data
772 if isfield(state, 'njobs')
773 self.njobs = state.njobs;
774 end
775 if isfield(state, 'ptaskcallers')
776 self.ptaskcallers = state.ptaskcallers;
777 end
778 if isfield(state, 'ilscaling')
779 self.ilscaling = state.ilscaling;
780 end
781
782 % Update layer models with imported state
783 it = 1;
784 if ~isempty(self.results)
785 it = size(self.results, 1);
786 end
787 self.updateLayers(it);
788
789 % Refresh all layer solvers with new parameters
790 for e = 1:self.nlayers
791 self.ensemble{e}.refreshChains();
792 switch self.solvers{e}.name
793 case {'SolverMVA', 'SolverNC'}
794 self.ensemble{e}.refreshRates();
795 otherwise
796 self.ensemble{e}.refreshProcesses();
797 end
798 self.solvers{e}.reset();
799 end
800 end
801
802 function update_solver(self, solverFactory)
803 % UPDATE_SOLVER Change the solver for all layers
804 %
805 % UPDATE_SOLVER(FACTORY) replaces all layer solvers with
806 % new solvers created by the given factory function.
807 %
808 % This allows switching between different solving methods
809 % (e.g., from MVA to JMT/DES) while preserving the current
810 % solution state.
811 %
812 % Example:
813 % solver = SolverLN(model, @(m) SolverMVA(m));
814 % solver.getEnsembleAvg(); % Fast initial solution
815 %
816 % solver.update_solver(@(m) SolverJMT(m, 'samples', 1e5));
817 % solver.getEnsembleAvg(); % Refine with simulation
818
819 self.solverFactory = solverFactory;
820
821 % Replace all layer solvers
822 for e = 1:self.nlayers
823 self.setSolver(solverFactory(self.ensemble{e}), e);
824 end
825 end
826
827 function [allMethods] = listValidMethods(self)
828 sn = self.model.getStruct();
829 % allMethods = LISTVALIDMETHODS()
830 % List valid methods for this solver
831 allMethods = {'default','moment3','mol'};
832 end
833 end
834
835 methods (Static)
836 function options = defaultOptions()
837 % OPTIONS = DEFAULTOPTIONS()
838 options = SolverOptions('LN');
839 end
840
841 function libs = getLibrariesUsed(sn, options)
842 % GETLIBRARIESUSED Get list of external libraries used by LN solver
843 % LN uses internal algorithms, no external library attribution needed
844 libs = {};
845 end
846 end
847end
848
849function solver = adaptiveSolverFactory(model, parentOptions)
850 % ADAPTIVESOLVERFACTORY - Select appropriate solver based on model characteristics
851 % Use JMT for models with open classes, MVA for pure closed networks
852 if nargin < 2
853 verbose = false;
854 else
855 verbose = parentOptions.verbose;
856 end
857
858 solver = SolverMVA(model, 'verbose', verbose);
859end