LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
refreshStruct.m
1function refreshStruct(self, hardRefresh)
2% REFRESHSTRUCT()
3%
4% Copyright (c) 2012-2026, Imperial College London
5% All rights reserved.
6
7sanitize(self);
8resolveSignals(self); % Resolve Signal placeholders to OpenSignal or ClosedSignal
9if nargin<2
10 hardRefresh = true;
11end
12
13%% store invariant information
14if self.hasStruct && ~hardRefresh
15 rtorig = self.sn.rtorig; % this must be destroyed with resetNetwork
16end
17
18if self.hasStruct && ~hardRefresh
19 nodetypes = sn.nodetypes;
20 classnames = sn.classnames;
21 nodenames = sn.nodenames;
22 refstat = sn.refstat;
23else
24 nodetypes = getNodeTypes(self);
25 classnames = getClassNames(self);
26 nodenames = getNodeNames(self);
27 refstat = getReferenceStations(self);
28
29 % Append FCR names and types to node lists (only when refreshing)
30 for f = 1:length(self.regions)
31 fcr = self.regions{f};
32 nodenames{end+1} = fcr.getName();
33 nodetypes(end+1) = NodeType.Region;
34 end
35end
36conn = self.getConnectionMatrix;
37njobs = getNumberOfJobs(self);
38numservers = getStationServers(self);
39lldscaling = getLimitedLoadDependence(self);
40cdscaling = getLimitedClassDependence(self);
41[ljdscaling, ljdcutoffs] = getLimitedJointDependence(self);
42[ljcdscaling, ljcdcutoffs] = getLimitedJointClassDependence(self);
43
44%% init minimal structure
45sn = NetworkStruct(); % create in self to ensure propagation
46if isempty(self.sn)
47 sn.rtorig = {};
48 sn.reward = {};
49else
50 sn.rtorig = self.sn.rtorig;
51 if isfield(self.sn, 'reward')
52 sn.reward = self.sn.reward; % preserve reward definitions
53 else
54 sn.reward = {};
55 end
56end
57% sn.nnodes counts physical nodes only; FCRs are virtual nodes appended to nodenames/nodetypes
58sn.nnodes = numel(self.nodes);
59sn.nclasses = length(classnames);
60
61%% get routing strategies
62routing = zeros(sn.nnodes, sn.nclasses);
63for ind=1:sn.nnodes
64 for r=1:sn.nclasses
65 if isempty(self.nodes{ind}.output.outputStrategy{r})
66 routing(ind,r) = RoutingStrategy.DISABLED;
67 else
68 routing(ind,r) = RoutingStrategy.fromText(self.nodes{ind}.output.outputStrategy{r}{2});
69 end
70 end
71end
72sn.isslc = false(sn.nclasses,1);
73for r=1:sn.nclasses
74 if isa(self.classes{r},'SelfLoopingClass')
75 sn.isslc(r) = true;
76 end
77end
78sn.issignal = false(sn.nclasses,1);
79sn.signaltype = cell(sn.nclasses,1);
80for r=1:sn.nclasses
81 sn.signaltype{r} = NaN;
82end
83% Detect Signal classes and populate issignal/signaltype
84% Check for Signal, OpenSignal, and ClosedSignal classes
85for r=1:sn.nclasses
86 if isa(self.classes{r},'Signal') || isa(self.classes{r},'OpenSignal') || isa(self.classes{r},'ClosedSignal')
87 sn.issignal(r) = true;
88 sn.signaltype{r} = self.classes{r}.signalType;
89 end
90end
91% Initialize syncreply - maps each class to its expected reply signal class index (-1 if none)
92sn.syncreply = -ones(sn.nclasses, 1);
93for r=1:sn.nclasses
94 if ~isempty(self.classes{r}.replySignalClass)
95 sn.syncreply(r) = self.classes{r}.replySignalClass.index - 1; % 0-based for JAR
96 end
97end
98% Initialize signal removal configuration fields
99sn.signalRemovalDist = cell(sn.nclasses,1);
100sn.signalRemovalPolicy = zeros(sn.nclasses, 1);
101sn.isCatastrophe = false(sn.nclasses, 1);
102% Populate removal configuration from Signal and ClosedSignal classes
103for r=1:sn.nclasses
104 if isa(self.classes{r}, 'Signal') || isa(self.classes{r}, 'OpenSignal')
105 if self.classes{r}.isCatastrophe()
106 sn.isCatastrophe(r) = true;
107 end
108 sn.signalRemovalDist{r} = self.classes{r}.removalDistribution;
109 sn.signalRemovalPolicy(r) = self.classes{r}.removalPolicy;
110 elseif isa(self.classes{r}, 'ClosedSignal')
111 sn.signalRemovalDist{r} = self.classes{r}.removalDistribution;
112 sn.signalRemovalPolicy(r) = self.classes{r}.removalPolicy;
113 end
114end
115sn.nclosedjobs = sum(njobs(isfinite(njobs)));
116sn.nservers = numservers;
117sn.isstation = (nodetypes == NodeType.Source | nodetypes == NodeType.Delay | nodetypes == NodeType.Queue | nodetypes == NodeType.Join | nodetypes == NodeType.Place);
118sn.nstations = sum(sn.isstation);
119sn.scv = ones(sn.nstations,sn.nclasses);
120sn.njobs = njobs(:)';
121sn.refstat = refstat;
122sn.space = cell(sn.nstations,1);
123sn.routing = routing;
124sn.chains = [];
125sn.lst = {};
126sn.lldscaling = lldscaling;
127sn.cdscaling = cdscaling;
128sn.ljdscaling = ljdscaling;
129sn.ljdcutoffs = ljdcutoffs;
130sn.ljcdscaling = ljcdscaling;
131sn.ljcdcutoffs = ljcdcutoffs;
132sn.nodetype = nodetypes;
133sn.nstations = sum(sn.isstation);
134sn.isstateful = (nodetypes == NodeType.Source | nodetypes == NodeType.Delay | nodetypes == NodeType.Queue | nodetypes == NodeType.Cache | nodetypes == NodeType.Join | nodetypes == NodeType.Router | nodetypes == NodeType.Place | nodetypes == NodeType.Transition);
135sn.isstatedep = false(sn.nnodes,3); % col 1: buffer, col 2: srv, col 3: routing
136sn.isfunction = [];
137for ind = 1:sn.nstations
138 if isa(self.stations{ind},'Queue')
139 sn.isfunction(ind) = ~isempty(self.stations{ind}.setupTime);
140 end
141end
142
143% Populate heterogeneous server fields from Queue nodes
144sn.nservertypes = zeros(sn.nstations, 1);
145sn.servertypenames = cell(sn.nstations, 1);
146sn.serverspertype = cell(sn.nstations, 1);
147sn.servercompat = cell(sn.nstations, 1);
148sn.heteroschedpolicy = zeros(sn.nstations, 1);
149for ist = 1:sn.nstations
150 if isa(self.stations{ist}, 'Queue') && ~isempty(self.stations{ist}.serverTypes)
151 nTypes = length(self.stations{ist}.serverTypes);
152 sn.nservertypes(ist) = nTypes;
153 sn.servertypenames{ist} = cell(1, nTypes);
154 sn.serverspertype{ist} = zeros(1, nTypes);
155 sn.servercompat{ist} = zeros(nTypes, sn.nclasses);
156
157 for t = 1:nTypes
158 st = self.stations{ist}.serverTypes{t};
159 sn.servertypenames{ist}{t} = st.getName();
160 sn.serverspertype{ist}(t) = st.numOfServers;
161
162 % Build compatibility matrix
163 for r = 1:sn.nclasses
164 if st.isCompatible(self.classes{r})
165 sn.servercompat{ist}(t, r) = 1;
166 end
167 end
168 end
169
170 % Get heterogeneous scheduling policy
171 if ~isempty(self.stations{ist}.heteroSchedPolicy)
172 sn.heteroschedpolicy(ist) = self.stations{ist}.heteroSchedPolicy;
173 end
174 end
175end
176
177for ind=1:sn.nnodes
178 switch sn.nodetype(ind)
179 case NodeType.Cache
180 sn.isstatedep(ind,2) = true; % state dependent service
181 % case NodeType.Place
182 % self.nodes{ind}.init();
183 % case NodeType.Transition
184 % self.nodes{ind}.init(); % this erases enablingConditions
185 end
186
187 for r=1:sn.nclasses
188 switch sn.routing(ind,r)
189 case {RoutingStrategy.RROBIN, RoutingStrategy.WRROBIN, RoutingStrategy.JSQ, RoutingStrategy.RL, RoutingStrategy.KCHOICES}
190 sn.isstatedep(ind,3) = true; % state dependent routing
191 end
192 end
193end
194
195sn.nstateful = sum(sn.isstateful);
196sn.state = cell(sn.nstations,1);
197for i=1:sn.nstateful
198 sn.state{i} = [];
199end
200sn.nodenames = nodenames;
201sn.classnames = classnames;
202sn.connmatrix = conn;
203
204sn.nodeToStateful =[];
205sn.nodeToStation =[];
206sn.stationToNode =[];
207sn.stationToStateful =[];
208sn.statefulToNode =[];
209sn.statefulToStation =[];
210for ind=1:sn.nnodes
211 sn.nodeToStateful(ind) = nd2sf(sn,ind);
212 sn.nodeToStation(ind) = nd2st(sn,ind);
213end
214for ist=1:sn.nstations
215 sn.stationToNode(ist) = st2nd(sn,ist);
216 sn.stationToStateful(ist) = st2sf(sn,ist);
217end
218for isf=1:sn.nstateful
219 sn.statefulToNode(isf) = sf2nd(sn,isf);
220 sn.statefulToStation(isf) = sf2st(sn,isf);
221end
222
223% Populate immediate feedback matrix (station x class)
224sn.immfeed = false(sn.nstations, sn.nclasses);
225for ist=1:sn.nstations
226 nodeIdx = sn.stationToNode(ist);
227 node = self.nodes{nodeIdx};
228 for r=1:sn.nclasses
229 % Check station-level setting (Queue only)
230 stationHas = false;
231 if isa(node, 'Queue') && ~isempty(node.immediateFeedback)
232 if ischar(node.immediateFeedback) && strcmp(node.immediateFeedback, 'all')
233 stationHas = true;
234 elseif iscell(node.immediateFeedback)
235 stationHas = any(cellfun(@(x) x == r, node.immediateFeedback));
236 end
237 end
238 % Check class-level setting
239 classHas = self.classes{r}.immediateFeedback;
240 sn.immfeed(ist, r) = stationHas || classHas;
241 end
242end
243
244sn.fj = self.getForkJoins();
245self.sn = sn;
246refreshPriorities(self);
247if exist('refreshDeadlines', 'file')
248 refreshDeadlines(self);
249else
250 % Inline implementation if method not loaded
251 K = getNumberOfClasses(self);
252 classdeadline = zeros(1,K);
253 for r=1:K
254 classdeadline(r) = self.getClassByIndex(r).deadline;
255 end
256 if ~isempty(self.sn)
257 self.sn.classdeadline = classdeadline;
258 end
259end
260refreshProcesses(self);
261
262% Export patience/impatience fields for abandonment-aware solvers (MAPMSG, JMT, etc.)
263sn = self.sn;
264sn.patienceProc = cell(sn.nstations, sn.nclasses);
265sn.impatienceClass = zeros(sn.nstations, sn.nclasses); % ImpatienceType (RENEGING, BALKING)
266sn.impatienceType = zeros(sn.nstations, sn.nclasses); % ProcessType (EXP, ERLANG, etc.)
267sn.impatienceMu = zeros(sn.nstations, sn.nclasses); % Rate parameter (1/mean)
268sn.impatiencePhi = zeros(sn.nstations, sn.nclasses); % SCV parameter
269sn.impatiencePhases = zeros(sn.nstations, sn.nclasses);
270sn.impatienceProc = cell(sn.nstations, sn.nclasses);
271sn.impatiencePie = cell(sn.nstations, sn.nclasses);
272for ist = 1:sn.nstations
273 node = self.stations{ist};
274 if isa(node, 'Queue')
275 for r = 1:sn.nclasses
276 patienceDist = node.getPatience(self.classes{r});
277 if ~isempty(patienceDist) && ~isa(patienceDist, 'Disabled')
278 % Convert patience distribution to MAP representation
279 patienceMAP = patienceDist.getProcess();
280 if iscell(patienceMAP) && length(patienceMAP) >= 2
281 sn.patienceProc{ist, r} = patienceMAP;
282 sn.impatienceProc{ist, r} = patienceMAP;
283
284 % Get the ProcessType of the patience distribution
285 if isprop(patienceDist, 'type')
286 sn.impatienceType(ist, r) = patienceDist.type;
287 elseif isa(patienceDist, 'Exp')
288 sn.impatienceType(ist, r) = ProcessType.EXP;
289 elseif isa(patienceDist, 'Erlang')
290 sn.impatienceType(ist, r) = ProcessType.ERLANG;
291 elseif isa(patienceDist, 'HyperExp')
292 sn.impatienceType(ist, r) = ProcessType.HYPEREXP;
293 elseif isa(patienceDist, 'Det')
294 sn.impatienceType(ist, r) = ProcessType.DET;
295 elseif isa(patienceDist, 'Gamma')
296 sn.impatienceType(ist, r) = ProcessType.GAMMA;
297 elseif isa(patienceDist, 'Pareto')
298 sn.impatienceType(ist, r) = ProcessType.PARETO;
299 elseif isa(patienceDist, 'Weibull')
300 sn.impatienceType(ist, r) = ProcessType.WEIBULL;
301 elseif isa(patienceDist, 'Lognormal')
302 sn.impatienceType(ist, r) = ProcessType.LOGNORMAL;
303 elseif isa(patienceDist, 'Uniform')
304 sn.impatienceType(ist, r) = ProcessType.UNIFORM;
305 elseif isa(patienceDist, 'Coxian')
306 sn.impatienceType(ist, r) = ProcessType.COXIAN;
307 elseif isa(patienceDist, 'APH')
308 sn.impatienceType(ist, r) = ProcessType.APH;
309 elseif isa(patienceDist, 'PH')
310 sn.impatienceType(ist, r) = ProcessType.PH;
311 else
312 % Default to PH for other Markovian distributions
313 sn.impatienceType(ist, r) = ProcessType.PH;
314 end
315
316 % Extract distribution parameters for JMT export
317 % Get rate (mu = 1/mean)
318 if ismethod(patienceDist, 'getMean')
319 meanVal = patienceDist.getMean();
320 if meanVal > 0
321 sn.impatienceMu(ist, r) = 1 / meanVal;
322 else
323 sn.impatienceMu(ist, r) = Inf;
324 end
325 else
326 % Estimate from MAP: mean = -pi * D0^(-1) * e
327 D0 = patienceMAP{1};
328 D1 = patienceMAP{2};
329 n = size(D0, 1);
330 pi = ones(1, n) / n; % Approximate stationary distribution
331 meanVal = -pi * (D0 \ ones(n, 1));
332 sn.impatienceMu(ist, r) = 1 / meanVal;
333 end
334
335 % Get SCV (phi)
336 if ismethod(patienceDist, 'getSCV')
337 sn.impatiencePhi(ist, r) = patienceDist.getSCV();
338 else
339 sn.impatiencePhi(ist, r) = 1.0; % Default to exponential SCV
340 end
341
342 % Get number of phases
343 if ismethod(patienceDist, 'getNumberOfPhases')
344 sn.impatiencePhases(ist, r) = patienceDist.getNumberOfPhases();
345 else
346 sn.impatiencePhases(ist, r) = size(patienceMAP{1}, 1);
347 end
348
349 % Get initial probability vector (pie)
350 n = size(patienceMAP{1}, 1);
351 D0 = patienceMAP{1};
352 D1 = patienceMAP{2};
353 % For PH: pie is from the MAP representation D1 = (-D0*e)*pie
354 exitRates = -D0 * ones(n, 1);
355 idx = find(exitRates > 1e-10, 1);
356 if ~isempty(idx)
357 sn.impatiencePie{ist, r} = D1(idx, :) / exitRates(idx);
358 else
359 sn.impatiencePie{ist, r} = ones(1, n) / n;
360 end
361 end
362 % Get impatience class (RENEGING, BALKING)
363 impType = node.getImpatienceType(self.classes{r});
364 if ~isempty(impType)
365 sn.impatienceClass(ist, r) = impType;
366 end
367 end
368 end
369 end
370end
371self.sn = sn;
372
373% Check if priorities are specified but no priority-aware scheduling policy is used
374sn = self.sn;
375if ~all(sn.classprio == sn.classprio(1))
376 % Priority classes exist, check if any station uses priority-aware scheduling
377 prioScheds = [SchedStrategy.PSPRIO, SchedStrategy.DPSPRIO, SchedStrategy.GPSPRIO, ...
378 SchedStrategy.HOL, SchedStrategy.FCFSPRIO, SchedStrategy.LCFSPRIO, ...
379 SchedStrategy.LCFSPRPRIO, SchedStrategy.LCFSPIPRIO];
380 if ~any(ismember(sn.sched, prioScheds))
381 line_warning(mfilename, 'Priority classes are specified but no priority-aware scheduling policy (PSPRIO, DPSPRIO, GPSPRIO, HOL, FCFSPRIO, LCFSPRIO, LCFSPRPRIO, LCFSPIPRIO) is used in the model. Priorities will be ignored.');
382 else
383 % Display priority info unless silent
384 global LINEVerbose;
385 if isempty(LINEVerbose) || LINEVerbose ~= VerboseLevel.SILENT
386 [minPrio, minIdx] = min(sn.classprio);
387 [maxPrio, maxIdx] = max(sn.classprio);
388 highestPrioClasses = find(sn.classprio == minPrio);
389 lowestPrioClasses = find(sn.classprio == maxPrio);
390 highNames = strjoin(arrayfun(@(i) sn.classnames{i}, highestPrioClasses, 'UniformOutput', false), ',');
391 lowNames = strjoin(arrayfun(@(i) sn.classnames{i}, lowestPrioClasses, 'UniformOutput', false), ',');
392 line_printf('Priority: highest=%s, lowest=%s\n', highNames, lowNames);
393 end
394 end
395end
396
397if any(nodetypes == NodeType.Cache)
398 % this also refreshes the routing matrix and the visits
399 refreshChains(self, false); % wantVisits
400else
401 % this also refreshes the routing matrix and the visits
402 refreshChains(self, true); % wantVisits
403end
404sn = self.sn;
405refclasses = getReferenceClasses(self);
406refclass = zeros(1,sn.nchains);
407for c=1:sn.nchains
408 isect = intersect(sn.inchain{c},find(refclasses));
409 if any(isect)
410 refclass(c) = isect;
411 end
412end
413sn.refclass = refclass;
414self.sn = sn;
415refreshLocalVars(self); % depends on chains (rtnodes)
416refreshPetriNetNodes(self);
417refreshSync(self); % this assumes that refreshChain is called before
418refreshGlobalSync(self);
419
420sn = self.sn;
421self.hasStruct = true;
422
423if any(sn.fj(:)) % if there are forks
424 % try to recompute visits after mixed-model transformation:
425 % we obtain the visits of auxiliary class from the transformed model
426 % then we sum them to the original classes
427 %try
428 [nonfjmodel, fjclassmap, forkmap, fanOut] = ModelAdapter.mmt(self);
429 if any(fanOut==1)
430 % in this case, one of the forks is degenerate with a single
431 % outgoing link so no visit correction is needed
432 line_warning(mfilename,'The specified fork-join topology has partial support, only SolverJMT simulation results may be reliable.\n');
433 % return;
434 end
435 fsn = nonfjmodel.getStruct();
436 for new_chain=(sn.nchains+1):fsn.nchains
437 anyAuxClass = fsn.inchain{new_chain}(1);
438 origFork = forkmap(anyAuxClass);
439 origChain = find(sn.chains(:,fjclassmap(anyAuxClass))); % original chain of the class
440 fsn.nodevisits{new_chain}(fsn.nodetype == NodeType.Source | fsn.nodetype == NodeType.Sink | fsn.nodetype == NodeType.Fork,:) = 0;
441 Vaux = fsn.nodevisits{new_chain}(:,fsn.inchain{new_chain});
442 if fsn.nnodes ~= sn.nnodes
443 % Build mapping from fsn nodes to sn nodes by name
444 % This handles cases where ClassSwitch/Source/Sink differ between models
445 VauxMapped = zeros(sn.nnodes, size(Vaux,2));
446 for fsnRow = 1:fsn.nnodes
447 nodeName = fsn.nodenames{fsnRow};
448 snRow = find(strcmp(sn.nodenames, nodeName), 1);
449 if ~isempty(snRow)
450 % This fsn node exists in sn - copy its visit data
451 VauxMapped(snRow, :) = Vaux(fsnRow, :);
452 end
453 % Nodes not in sn (ClassSwitch, Source, Sink) are skipped
454 end
455 Vaux = VauxMapped;
456 end
457 X = sn.nodevisits{origChain};
458 for jaux=1:length(fsn.inchain{new_chain})
459 j = sn.inchain{origChain}(jaux); % class in the original chain
460 self.sn.nodevisits{origChain}(:,j) = sn.nodeparam{origFork}.fanOut*(X(:,j) + Vaux(:,jaux));
461 end
462 end
463end
464
465sn = refreshRegions(self);
466self.sn = sn;
467end
468
469function stat_idx = nd2st(sn, node_idx)
470% STAT_IDX = ND2ST(NODE_IDX)
471
472if sn.isstation(node_idx)
473 stat_idx = at(cumsum(sn.isstation),node_idx);
474else
475 stat_idx = NaN;
476end
477end
478
479function node_idx = st2nd(sn,stat_idx)
480% NODE_IDX = ST2ND(SELF,STAT_IDX)
481
482v = cumsum(sn.isstation) == stat_idx;
483if any(v)
484 node_idx = find(v, 1);
485else
486 node_idx = NaN;
487end
488end
489
490function sful_idx = st2sf(sn,stat_idx)
491% SFUL_IDX = ST2SF(SELF,STAT_IDX)
492
493sful_idx = nd2sf(sn,st2nd(sn,stat_idx));
494end
495
496function sful_idx = nd2sf(sn, node_idx)
497% SFUL_IDX = ND2SF(NODE_IDX)
498
499if sn.isstateful(node_idx)
500 sful_idx = at(cumsum(sn.isstateful),node_idx);
501else
502 sful_idx = NaN;
503end
504end
505
506function node_idx = sf2nd(sn,stat_idx)
507% NODE_IDX = SF2ND(SELF,STAT_IDX)
508
509v = cumsum(sn.isstateful) == stat_idx;
510if any(v)
511 node_idx = find(v, 1);
512else
513 node_idx = NaN;
514end
515end
516
517function stat_idx = sf2st(sn,sful_idx)
518% STAT_IDX = SF2ST(SELF,SFUL_IDX)
519
520stat_idx = nd2st(sn,sf2nd(sn,sful_idx));
521end
Definition mmt.m:92