LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
solver_ssa_nrm.m
1function [QN, UN, RN, TN, CN, XN, lG, sn] = solver_ssa_nrm(sn, options)
2% SOLVER_SSA_NRM Steady‑state analysis via the Next‑Reaction Method (SSA)
3%
4% [QN, UN, RN, TN, CN, XN, LG, SN] = SOLVER_SSA_NRM(SN, OPTIONS)
5% runs a stochastic simulation of the queueing network described in SN
6% for OPTIONS.samples reaction firings using Gibson & Bruck's
7% Next‑Reaction Method. During the run it:
8% • computes performance metrics directly during simulation;
9% • returns standard queueing performance measures.
10%
11% Outputs
12% QN – M×K matrix of mean queue lengths
13% UN – M×K matrix of utilizations
14% RN – M×K matrix of response times
15% TN – M×K matrix of throughputs
16% CN – 1×K vector of cycle times
17% XN – 1×K vector of system throughputs
18% LG – Logarithm of normalizing constant (not computed)
19% SN – (Possibly updated) network structure.
20%
21% See also SOLVER_SSA_NRM_SPACE, NEXT_REACTION_METHOD_DIRECT.
22
23% ---------------------------------------------------------------------
24% Parameters & shorthands
25% ---------------------------------------------------------------------
26samples = options.samples;
27R = sn.nclasses;
28I = sn.nnodes;
29M = sn.nstations;
30K = sn.nclasses;
31state = sn.state;
32
33% ---------------------------------------------------------------------
34% Stoichiometry & reaction mapping (self‑loops included) ----------------
35% ---------------------------------------------------------------------
36S = zeros(0, I*R); % will transpose at the end
37fromIdx = [];
38toIdx = [];
39fromIR = [];
40
41% this is currently M^2*R^2, it can be lowered to M*R decoupling the
42% routing
43k = 0;
44for ind = 1:I
45 for r = 1:R
46 k = k + 1;
47 fromIR(k,:) = [ind, r];
48 fromIdx(k) = (ind-1)*R + r;
49 probIR{k} = [];
50 toIdx{k} = [];
51 Srow = zeros(1, I*R); % build stoichiometry row
52 if sn.isslc(r)
53 Srow(fromIdx(k)) = -Inf;
54 else
55 Srow(fromIdx(k)) = -1;
56 for jnd = 1:I
57 for s = 1:R
58 if sn.rtnodes((ind-1)*R+r, (jnd-1)*R+s) > 0
59 toIdx{k}(end+1) = (jnd-1)*R + s;
60 p = sn.rtnodes((ind-1)*R+r, (jnd-1)*R+s);
61 probIR{k}(end+1) = p;
62 Srow((jnd-1)*R + s) = Srow((jnd-1)*R + s) + p;
63 end
64 end
65 end
66 end
67 S(k,:) = Srow;
68 end
69end
70S = S.'; % states × reactions
71
72% ---------------------------------------------------------------------
73% Initial state vector --------------------------------------------------
74% ---------------------------------------------------------------------
75nvec0 = zeros(I*R,1); % initial state (aggregate state)
76for ind=1:I
77 if sn.isstateful(ind)
78 [~,nir] = State.toMarginalAggr(sn, ind, state{sn.nodeToStateful(ind)});
79 for r = 1:R
80 if isinf(nir(r))
81 if sn.nodetype(ind) == NodeType.Source
82 nir(r) = 1;
83 else
84 line_error(mfilename, 'Infinite population error.');
85 end
86 end
87 nvec0((ind-1)*R + r,1) = nir(r);
88 end
89 end
90end
91
92mi = zeros(I,1);
93rates = zeros(I,R);
94for ind=1:I
95 if sn.isstation(ind)
96 for r=1:R
97 ist = sn.nodeToStation(ind);
98 muir = sn.rates(ist,r);
99 if ~isnan(muir)
100 rates(ind,r) = muir;
101 end
102 mi(ind,1) = sn.nservers(ist);
103 end
104 else
105 for r=1:R
106 rates(ind,r) = GlobalConstants.Immediate;
107 mi(ind,1) = GlobalConstants.MaxInt;
108 end
109 end
110 mi(isinf(mi)) = GlobalConstants.MaxInt;
111end
112
113% Propensity function ---------------------------------------------------
114epstol = GlobalConstants.Zero;
115a = {};
116for j=1:length(fromIdx)
117 if sn.isstation(fromIR(j,1))
118 switch sn.sched(sn.nodeToStation(fromIR(j,1)))
119 case SchedStrategy.EXT
120 a{j} = @(X) rates(fromIR(j,1), fromIR(j,2));
121 case SchedStrategy.INF
122 a{j} = @(X) rates(fromIR(j,1), fromIR(j,2)) * X(fromIdx(j));
123 case SchedStrategy.PS
124 if R == 1 % single class
125 a{j} = @(X) rates(fromIR(j,1), fromIR(j,2)) * min( mi(fromIR(j,1)), X(fromIdx(j)));
126 else
127 a{j} = @(X) rates(fromIR(j,1), fromIR(j,2)) * ( X(fromIdx(j)) ./ ...
128 (epstol+sum( X(((fromIR(j,1)-1)*R + 1):((fromIR(j,1)-1)*R + R)) ) )) * ...
129 min( mi(fromIR(j,1)), ...
130 (epstol+sum( X(((fromIR(j,1)-1)*R + 1):((fromIR(j,1)-1)*R + R)) ) ));
131 end
132 end
133 else
134 a{j} = @(X) rates(fromIR(j,1), fromIR(j,2)) * min(1, X(fromIdx(j)));
135 end
136end
137
138% Propensity functions dependencies -----------------------------------
139D = cell(1,size(S,2));
140for k=1:size(D,2)
141 J = find(S(:,k))'; % set of state variables affected by reaction k
142 vecd = [];
143 for j=1:length(J)
144 % (ind-1)*R + r
145 pos = J(j);
146 r = mod(pos-1, R) + 1;
147 ind = ((pos-r)/R) + 1;
148 vecd(end+1:end+R) = ((ind-1)*R + 1) : (ind*R);
149 end
150 % vecd now contains all state variables affected by the firing of
151 % reaction k. We now find the propensity functions that depend
152 % on those variables
153 if ~isempty(vecd)
154 vecd = unique(vecd);
155 vecs = [];
156 for j=1:length(vecd)
157 vecs = [vecs,find(S(vecd(j),:)<0)];
158 end
159 D{k} = unique(vecs);
160 else
161 D{k} = [];
162 end
163end
164
165% Having accounted for them in D, we can now remove self-loops markings
166S(isinf(S))=0;
167
168% ---------------------------------------------------------------------
169% Initialize performance metric matrices
170% ---------------------------------------------------------------------
171lG = 0; % Not computed in SSA
172
173% ---------------------------------------------------------------------
174% Run SSA/NRM with direct metric computation
175% ---------------------------------------------------------------------
176[QN, UN, RN, TN, CN, XN] = next_reaction_method_direct(S, D, a, nvec0, samples, options, sn, fromIdx, fromIR);
177
178end % solver_ssa_nrm
179
180% ======================================================================
181% Next-Reaction Method with direct metric computation
182% ======================================================================
183function [QN, UN, RN, TN, CN, XN] = next_reaction_method_direct(S, D, a, nvec0, samples, options, sn, fromIdx, fromIR)
184
185numReactions = size(S,2);
186R = sn.nclasses;
187I = sn.nnodes;
188M = sn.nstations;
189K = sn.nclasses;
190QN = zeros(M, K);
191UN = zeros(M, K);
192RN = zeros(M, K);
193TN = zeros(M, K);
194CN = zeros(1, K);
195XN = zeros(1, K);
196
197% when a reaction fires, this matrix helps selecting the probability that a
198% particular routing or phase is selected as a result ------------------
199P = S; P(P<0)=P(P<0)+1';
200fromIdxCell = cell(numReactions,1);
201toIdxCell = cell(numReactions,1);
202cdfVec = cell(numReactions,1);
203for r=1:numReactions
204 nnzP(r) = nnz(P(:,r));
205 if nnzP(r)>1
206 fromIdxCell{r} = find(S(:,r)<0);
207 toIdxCell{r} = find(P(:,r));
208 cdfVec{r} = cumsum(P(toIdxCell{r},r));
209 end
210end
211
212% initialise Gillespie clocks ------------------------------------------
213t = 0;
214for k=1:size(S,2)
215 Ak(k) = a{k}(nvec0);
216end
217nvec = nvec0;
218Pk = -log(rand(1,numReactions));
219Tk = zeros(1,numReactions);
220
221tau = (Pk - Tk) ./ Ak;
222
223% Performance tracking variables
224totalTime = 0;
225NK = sn.njobs'; % Jobs per class
226servers = sn.nservers;
227
228n = 1;
229while n <= samples
230 [dt, kfire] = min(tau);
231 if isinf(dt), line_error(mfilename,'Deadlock. Quitting nrm method.'); end
232
233 totalTime = totalTime + dt;
234
235 % Accumulate state-dependent metrics during this time interval
236 for ist = 1:M
237 ind = sn.stationToNode(ist);
238 for k = 1:K
239 currentPop = nvec((ind-1)*R + k);
240
241 % Accumulate queue length (QN)
242 QN(ist, k) = QN(ist, k) + currentPop * dt;
243
244 % Compute throughput contribution from departures
245 irIndex = (ind-1)*R + k;
246 depRate = 0;
247 idx = find(fromIdx == irIndex);
248 if ~isempty(idx)
249 depRate = Ak(idx(1));
250 end
251 TN(ist, k) = TN(ist, k) + depRate * dt;
252
253 % Compute utilization based on scheduling policy
254 switch sn.sched(ist)
255 case {SchedStrategy.INF, SchedStrategy.EXT}
256 UN(ist, k) = UN(ist, k) + currentPop * dt;
257 case SchedStrategy.PS
258 totalPop = sum(nvec(((ind-1)*R + 1):(ind*R)));
259 if totalPop > 0
260 utilization = (currentPop / totalPop) * min(servers(ist), totalPop) / servers(ist);
261 else
262 utilization = 0;
263 end
264 UN(ist, k) = UN(ist, k) + utilization * dt;
265 end
266 end
267 end
268
269 t = t + dt;
270
271 % update aggregate state
272 if nnzP(kfire)>1
273 r = 1+find(rand>=cdfVec{kfire},1);
274 if isempty(r)
275 r = 1;
276 end
277 nvec(fromIdxCell{kfire}) = nvec(fromIdxCell{kfire}) - 1;
278 nvec(toIdxCell{kfire}(r)) = nvec(toIdxCell{kfire}(r)) + 1;
279 else
280 nvec = nvec + S(:,kfire); % zero change for self-loops
281 end
282
283 Tk = Tk + Ak * dt;
284
285 % update rates for all reactions dependent on the last fired reaction
286 for k=D{kfire}
287 Ak(k) = a{k}(nvec);
288 end
289
290 % update clocks
291 Pk(kfire) = Pk(kfire) - log(rand);
292 tau = (Pk - Tk) ./ Ak;
293 tau(Ak==0) = inf;
294
295 % do not count immediate events
296 n = n + 1;
297 print_progress(options, n);
298end % while
299% Print newline after progress counter
300if isfield(options,'verbose') && options.verbose
301 line_printf('\n');
302end
303
304% Normalize metrics by total time
305if totalTime > 0
306 for ist = 1:M
307 for k = 1:K
308 QN(ist, k) = QN(ist, k) / totalTime;
309 UN(ist, k) = UN(ist, k) / totalTime;
310 TN(ist, k) = TN(ist, k) / totalTime;
311 end
312 end
313end
314
315% Compute derived metrics
316for k = 1:K
317 % System throughput at reference station
318 XN(1, k) = TN(sn.refstat(k), k);
319
320 % Response times
321 for ist = 1:M
322 if TN(ist, k) > 0
323 RN(ist, k) = QN(ist, k) / TN(ist, k);
324 else
325 RN(ist, k) = 0;
326 end
327 end
328
329 % Cycle times
330 if XN(1, k) > 0
331 CN(1, k) = NK(k) / XN(1, k);
332 end
333end
334
335% Handle NaN values
336QN(isnan(QN)) = 0;
337UN(isnan(UN)) = 0;
338RN(isnan(RN)) = 0;
339XN(isnan(XN)) = 0;
340TN(isnan(TN)) = 0;
341CN(isnan(CN)) = 0;
342
343 function print_progress(opt, samples_collected)
344 if ~isfield(opt,'verbose') || ~opt.verbose, return; end
345 if samples_collected == 1e3
346 line_printf('\nSSA samples: %8d', samples_collected);
347 elseif opt.verbose == 2
348 if samples_collected == 0
349 line_printf('\nSSA samples: %9d', samples_collected);
350 else
351 line_printf('\b\b\b\b\b\b\b\b\b%9d', samples_collected);
352 end
353 elseif mod(samples_collected,1e3)==0 || opt.verbose == 2
354 line_printf('\b\b\b\b\b\b\b\b\b%9d', samples_collected);
355 end
356 end
357end % next_reaction_method_direct