LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
solver_fluid_matrix.m
1function [QN,UN,RN,TN,xvec_it,QNt,UNt,TNt,xvec_t,t,iters,runtime] = solver_fluid_matrix(sn, options)
2
3% [QN,UN,RN,TN,CN,RUNTIME] = SOLVER_FLUID_MATRIX(QN, OPTIONS)
4
5% Copyright (c) 2012-2026, Imperial College London
6% All rights reserved.
7
8M = sn.nstations; %number of stations
9K = sn.nclasses; %number of classes
10pie = sn.pie;
11PH = sn.proc;
12P_full = sn.rt; % Full routing matrix (stateful nodes)
13NK = sn.njobs'; %initial population
14S = sn.nservers;
15infServers = isinf(S);
16S(infServers) = sum(NK);
17nphases = sn.phases;
18%refstat = sn.refstat; % reference station
19weights = ones(M,K);
20
21% Extract station-to-station routing matrix from stateful-to-stateful matrix
22% using stochastic complementation to resolve routing through non-station
23% stateful nodes (e.g., Router nodes)
24station_indices = [];
25for ist = 1:M
26 isf = sn.stationToStateful(ist);
27 for r = 1:K
28 station_indices = [station_indices, (isf-1)*K + r];
29 end
30end
31P = dtmc_stochcomp(P_full, station_indices);
32
33% Remove Sink->Source feedback routing for open classes
34% In open networks, jobs exit at Sink and should not recirculate back to Source.
35% The routing matrix includes this feedback (added by getRoutingMatrix.m for
36% pseudo-closed network analysis), but it causes incorrect flow balance in the
37% fluid ODE formulation because arrivals are already accounted for via Alambda.
38for src_ist = 1:M
39 if sn.sched(src_ist) == SchedStrategy.EXT
40 % This is a Source station - remove feedback routing TO it for open classes
41 for r = 1:K
42 % Check if this is an open class with external arrivals
43 if ~isnan(sn.rates(src_ist, r)) && sn.rates(src_ist, r) > 0
44 % Zero out routing TO this Source for this class from all other stations
45 src_col = (src_ist - 1) * K + r; % Column index in P for (Source, class r)
46 for from_ist = 1:M
47 if from_ist ~= src_ist % Don't modify Source's own outgoing routing
48 for from_r = 1:K
49 from_row = (from_ist - 1) * K + from_r;
50 P(from_row, src_col) = 0; % Remove feedback to Source
51 end
52 end
53 end
54 end
55 end
56 end
57end
58
59% ODE building as per Ruuskanen et al., PEVA 151 (2021).
60Psi = [];
61A = [];
62B = [];
63for ist=1:M
64 for r=1:K
65 if nphases(ist,r)==0
66 Psi = blkdiag(Psi,0);
67 B = blkdiag(B,0);
68 A = blkdiag(A,NaN);
69 else
70 Psi = blkdiag(Psi,PH{ist}{r}{1});
71 B = blkdiag(B,sum(PH{ist}{r}{2},2));
72 A = blkdiag(A,pie{ist}{r}');
73 end
74 end
75end
76W = Psi + B*P*A';
77
78% Build arrival rate vector Aλ for mixed/open networks (Ruuskanen et al., PEVA 151 (2021), Eq. 7)
79% Following the ground truth implementation: Source is conceptually excluded from state space.
80% Arrivals go directly to queue phases (not Source phases) weighted by routing probabilities.
81% This matches dx = W' * theta + A * lambda where lambda represents arrivals INTO queues.
82
83% First, identify Source stations and their arrival rates per class
84source_arrivals = zeros(M, K); % source_arrivals(src, r) = arrival rate at source src for class r
85for src_ist = 1:M
86 if sn.sched(src_ist) == SchedStrategy.EXT
87 for r = 1:K
88 if ~isnan(sn.rates(src_ist, r)) && sn.rates(src_ist, r) > 0
89 source_arrivals(src_ist, r) = sn.rates(src_ist, r);
90 end
91 end
92 end
93end
94
95% Build Alambda: arrivals go to QUEUE phases (where jobs route from Source), not Source phases
96Alambda_full = zeros(size(A,1), 1);
97state = 0;
98for ist=1:M
99 for r=1:K
100 if nphases(ist,r) > 0
101 if sn.sched(ist) == SchedStrategy.EXT
102 % Source station: do NOT add arrivals here (arrivals go to downstream queues)
103 state = state + nphases(ist,r);
104 else
105 % Queue station: check if it receives arrivals from any Source
106 arrival_rate_to_queue = 0;
107 for src_ist = 1:M
108 if source_arrivals(src_ist, r) > 0
109 % Get routing probability from Source to this queue for this class
110 src_row = (src_ist - 1) * K + r;
111 queue_col = (ist - 1) * K + r;
112 routing_prob = P(src_row, queue_col);
113 arrival_rate_to_queue = arrival_rate_to_queue + source_arrivals(src_ist, r) * routing_prob;
114 end
115 end
116
117 if arrival_rate_to_queue > 0
118 % Apply arrivals according to entrance probability ζ^{i,r}
119 for k=1:nphases(ist,r)
120 state = state + 1;
121 Alambda_full(state) = pie{ist}{r}(k) * arrival_rate_to_queue;
122 end
123 else
124 state = state + nphases(ist,r);
125 end
126 end
127 else
128 % Add placeholder for disabled class to match W matrix structure
129 state = state + 1;
130 % Alambda_full(state) stays 0 since there are no arrivals
131 end
132 end
133end
134
135% remove disabled transitions
136keep = find(~isnan(sum(W,1)));
137W = W(keep,:);
138W = W(:,keep);
139Alambda = Alambda_full(keep); % Also filter arrival vector
140
141% Eliminate immediate transitions if requested
142state_map_imm = [];
143W_pre_elim = [];
144Alambda_pre_elim = [];
145if isfield(options.config, 'hide_immediate') && options.config.hide_immediate
146 W_pre_elim = W; % Save for Alambda correction
147 Alambda_pre_elim = Alambda;
148 [W, state_map_imm] = eliminate_immediate_matrix(W, sn, options);
149end
150
151Qa = []; % state mapping to queues (called Q(a) in Ruuskanen et al.)
152SQC = zeros(M*K,0); % to compute per-class queue length at the end
153SUC = zeros(M*K,0); % to compute per-class utilizations at the end
154STC = zeros(M*K,0); % to compute per-class throughput at the end
155x0_build = []; % Build x0 with same structure as W matrix
156%x0 = []; % initial state
157state = 0;
158init_sol_idx = 0; % Index into options.init_sol for enabled classes
159for ist=1:M
160 for r=1:K
161 if nphases(ist,r)==0
162 % Add placeholder for disabled transition (matching W matrix structure)
163 state = state + 1;
164 Qa(1,state) = ist; %#ok<*AGROW>
165 SQC((ist-1)*K+r,state) = 0; % No queue contribution
166 SUC((ist-1)*K+r,state) = 0; % No utilization contribution
167 STC((ist-1)*K+r,state) = 0; % No throughput contribution
168 x0_build(state,1) = 0; % No initial population
169 else
170 for k=1:nphases(ist,r)
171 state = state + 1;
172 init_sol_idx = init_sol_idx + 1;
173 Qa(1,state) = ist;
174 SQC((ist-1)*K+r,state) = 1;
175 SUC((ist-1)*K+r,state) = 1/S(ist);
176 STC((ist-1)*K+r,state) = sum(sn.proc{ist}{r}{2}(k,:));
177 x0_build(state,1) = options.init_sol(init_sol_idx);
178 end
179 end
180 % code to initialize all jobs at ref station
181 %if i == refstat(r)
182 % x0 = [x0; NK(r)*pie{i}{r}']; % initial probability of PH
183 %else
184 % x0 = [x0; zeros(nphases(i,r),1)];
185 %end
186 end
187end
188x0 = x0_build;
189
190% Apply keep filtering for disabled transitions (NaN in W matrix)
191Qa = Qa(keep);
192SQC = SQC(:, keep);
193SUC = SUC(:, keep);
194STC = STC(:, keep);
195x0 = x0(keep);
196
197% Save full matrices before immediate elimination for metric reconstruction
198Qa_full = Qa;
199STC_full = STC;
200imm_states_in_keep = []; % Indices of immediate states within keep-filtered space
201
202% Apply state mapping if immediate transitions were eliminated
203if ~isempty(state_map_imm)
204 % Identify eliminated (immediate) states
205 imm_states_in_keep = setdiff(1:length(Qa), state_map_imm);
206
207 Qa = Qa(state_map_imm);
208 SQC = SQC(:, state_map_imm);
209 SUC = SUC(:, state_map_imm);
210 STC = STC(:, state_map_imm);
211 x0 = x0(state_map_imm);
212
213 % Correct Alambda for arrivals directed at eliminated immediate states.
214 % From quasi-steady-state assumption (dx_I/dt = 0):
215 % lambda_red = lambda_T - Q_IT' * (Q_II')^{-1} * lambda_I
216 % where T = timed states, I = immediate states.
217 lambda_T = Alambda_pre_elim(state_map_imm);
218 lambda_I = Alambda_pre_elim(imm_states_in_keep);
219 if any(lambda_I ~= 0)
220 Q_IT = W_pre_elim(imm_states_in_keep, state_map_imm);
221 Q_II = W_pre_elim(imm_states_in_keep, imm_states_in_keep);
222 Alambda = lambda_T - Q_IT' * (Q_II' \ lambda_I);
223 else
224 Alambda = lambda_T;
225 end
226end
227
228% Build SQ matrix to compute total queue length per station in ODEs
229% SQ(s,:) sums all states at the same station as state s
230nstates = length(x0);
231SQ = zeros(nstates, nstates);
232for s = 1:nstates
233 ist = Qa(s);
234 SQ(s, Qa == ist) = 1; % Sum all states at the same station
235end
236
237% Identify Source station states (EXT scheduler)
238% For Source stations, theta should be 0.0 to effectively bypass Source in dynamics.
239% This matches the ground truth implementation where Source is excluded from state space.
240% Arrivals are injected directly into queue phases via Alambda.
241isSourceState = false(nstates, 1);
242for s = 1:nstates
243 ist = Qa(s);
244 if sn.sched(ist) == SchedStrategy.EXT
245 isSourceState(s) = true;
246 x0(s) = 0; % Initialize Source phases to 0 (no mass at Source)
247 end
248end
249
250%x0
251
252tol = options.tol;
253timespan = options.timespan;
254itermax = options.iter_max;
255odeopt = odeset('AbsTol', tol, 'RelTol', tol, 'NonNegative', 1:length(x0));
256nonZeroRates = abs(W(abs(W)>0)); nonZeroRates=nonZeroRates(:);
257trange = [timespan(1),min(timespan(2),abs(10*itermax/min(nonZeroRates)))];
258
259% Check if p-norm smoothing should be used (pstar parameter)
260use_pnorm = isfield(options, 'pstar') && ~isempty(options.pstar) || ...
261 (isfield(options.config, 'pstar') && ~isempty(options.config.pstar));
262
263if use_pnorm
264 % Get pstar values - expand scalar to per-station array
265 if isfield(options, 'pstar') && ~isempty(options.pstar)
266 pstar_val = options.pstar;
267 else
268 pstar_val = options.config.pstar;
269 end
270 if isscalar(pstar_val)
271 pstar_val = pstar_val * ones(M, 1);
272 end
273 % Create per-phase pstar array (pQa) using the filtered Qa mapping
274 pQa = pstar_val(Qa(:)); % Use Qa which is already filtered by keep and state_map_imm
275 Sa_pnorm = S(Qa(:)); % Column vector for pnorm_ode
276end
277
278T0 = tic;
279iters = 1;
280if use_pnorm
281 % p-norm smoothing ODE as per Ruuskanen et al., PEVA 151 (2021)
282 % ghat = 1 / (1 + (x/c)^p)^(1/p) where x is queue length, c is servers, p is pstar
283 % dx/dt = W^T * θ̂(x,p) + Aλ (Eq. 27 for mixed networks)
284 ode_pnorm_func = @(t,x) pnorm_ode(x, W, SQ, Sa_pnorm, pQa, Alambda, isSourceState);
285 if options.stiff
286 [t, xvec_t] = ode_solve_stiff(ode_pnorm_func, trange, x0, odeopt, options);
287 else
288 [t, xvec_t] = ode_solve(ode_pnorm_func, trange, x0, odeopt, options);
289 end
290else
291 % Standard matrix method without smoothing
292 % dx/dt = W^T * θ(x) + Aλ (Eq. 12 for mixed networks)
293 Sa_ode = S(Qa(:)); % Column vector for element-wise operations
294 % Define theta function with special handling for Source stations
295 theta_func = @(x) compute_theta(x, SQ, Sa_ode, isSourceState);
296 if options.stiff
297 [t, xvec_t] = ode_solve_stiff(@(t,x) W'*theta_func(x) + Alambda, trange, x0, odeopt, options);
298 else
299 [t, xvec_t] = ode_solve(@(t,x) W'*theta_func(x) + Alambda, trange, x0, odeopt, options);
300 end
301end
302runtime = toc(T0);
303
304Tmax = size(xvec_t,1);
305QNtmp = cell(1,Tmax);
306UNtmp = cell(1,Tmax);
307RNtmp = cell(1,Tmax);
308TNtmp = cell(1,Tmax);
309Sa = S(Qa(:)); % Column vector for element-wise operations
310S = repmat(S,1,K)'; S=S(:);
311for j=1:Tmax
312 x = xvec_t(j,:)';
313 QNtmp{j} = zeros(K,M);
314 TNtmp{j} = zeros(K,M);
315 UNtmp{j} = zeros(K,M);
316 RNtmp{j} = zeros(K,M);
317
318 QNtmp{j}(:) = SQC*x;
319 % Use compute_theta for consistent handling of Source stations
320 theta_j = compute_theta(x, SQ, Sa, isSourceState);
321 TNtmp{j}(:) = STC*theta_j;
322 UNtmp{j}(:) = SUC*theta_j;
323 % Little's law is invalid in transient so this vector is not returned
324 % except the last element as an approximation of the actual RN
325 RNtmp{j}(:) = QNtmp{j}(:)./TNtmp{j}(:);
326
327 QNtmp{j} = QNtmp{j}';
328 UNtmp{j} = UNtmp{j}';
329 RNtmp{j} = RNtmp{j}';
330 TNtmp{j} = TNtmp{j}';
331end
332% steady state metrics
333for j=1:Tmax
334 QNtmp{j} = QNtmp{j}(:);
335 UNtmp{j} = UNtmp{j}(:);
336 RNtmp{j} = RNtmp{j}(:);
337 TNtmp{j} = TNtmp{j}(:);
338end
339
340% compute cell array with time-varying metrics for stations and classes
341QNtmp = cell2mat(QNtmp)';
342UNtmp = cell2mat(UNtmp)';
343RNtmp = cell2mat(RNtmp)';
344TNtmp = cell2mat(TNtmp)';
345QNt = cell(M,K);
346UNt = cell(M,K);
347RNt = cell(M,K);
348TNt = cell(M,K);
349for ist=1:M
350 for r=1:K
351 QNt{ist,r} = QNtmp(:,(r-1)*M+ist);
352 UNt{ist,r} = UNtmp(:,(r-1)*M+ist);
353 RNt{ist,r} = RNtmp(:,(r-1)*M+ist);
354 TNt{ist,r} = TNtmp(:,(r-1)*M+ist);
355 end
356end
357QN = reshape(QNtmp(end,:),M,K);
358UN = reshape(UNtmp(end,:),M,K);
359RN = reshape(RNtmp(end,:),M,K);
360TN = reshape(TNtmp(end,:),M,K);
361
362% Set throughput at Source stations to arrival rates for open classes
363% Source stations have theta = 0 in the ODE (to bypass Source in dynamics),
364% but their throughput should equal the external arrival rate.
365for ist=1:M
366 if sn.sched(ist) == SchedStrategy.EXT
367 for r=1:K
368 if ~isnan(sn.rates(ist, r)) && sn.rates(ist, r) > 0
369 TN(ist, r) = sn.rates(ist, r);
370 end
371 end
372 end
373end
374
375% Reconstruct throughput for eliminated immediate states using flow conservation
376if ~isempty(imm_states_in_keep)
377 % For each eliminated state, throughput = sum of incoming throughputs via routing
378 for idx = 1:length(imm_states_in_keep)
379 s = imm_states_in_keep(idx);
380 ist = Qa_full(s); % Station of this eliminated state
381 % Find which class this state belongs to
382 r = 0;
383 state_count = 0;
384 for rr = 1:K
385 state_count = state_count + nphases(ist, rr);
386 if s <= sum(Qa_full == ist & (1:length(Qa_full)) <= state_count)
387 r = rr;
388 break;
389 end
390 end
391 if r == 0
392 % Fallback: find class by examining STC_full
393 for rr = 1:K
394 if STC_full((ist-1)*K+rr, s) > 0
395 r = rr;
396 break;
397 end
398 end
399 end
400 if r > 0
401 % Compute incoming throughput using routing matrix P
402 % T(ist, r) = sum over all (j, l) of T(j, l) * P((j,l) -> (ist,r))
403 incoming_tput = 0;
404 for j = 1:M
405 for l = 1:K
406 p_jl_ir = P((j-1)*K+l, (ist-1)*K+r);
407 if p_jl_ir > 0
408 incoming_tput = incoming_tput + TN(j, l) * p_jl_ir;
409 end
410 end
411 end
412 % Also add external arrivals if this is a source station
413 if sn.sched(ist) == SchedStrategy.EXT && ~isnan(sn.rates(ist, r))
414 incoming_tput = incoming_tput + sn.rates(ist, r);
415 end
416 TN(ist, r) = incoming_tput;
417 end
418 end
419end
420
421xvec_it = {xvec_t(end,:)};
422end
423
424function dxdt = pnorm_ode(x, W, SQ, Sa, pQa, Alambda, isSourceState)
425% PNORM_ODE - ODE derivative using p-norm smoothing
426% As per Ruuskanen et al., PEVA 151 (2021)
427% dxdt = W' * (x .* ghat) + Aλ (Eq. 27)
428% where ghat = 1 / (1 + (sumXQa/Sa)^pQa)^(1/pQa)
429
430sumXQa = GlobalConstants.FineTol + SQ * x;
431ghat = zeros(size(x));
432for i = 1:length(x)
433 xVal = sumXQa(i);
434 cVal = Sa(i);
435 pVal = pQa(i);
436 if cVal > 0 && pVal > 0
437 ghatVal = 1.0 / (1 + (xVal / cVal)^pVal)^(1/pVal);
438 if isnan(ghatVal)
439 ghat(i) = 0;
440 else
441 ghat(i) = ghatVal;
442 end
443 else
444 ghat(i) = 1;
445 end
446end
447
448% Compute effective rate (x .* ghat), with special handling for Source stations
449theta_eff = x .* ghat;
450% For Source stations, override to 0.0 to bypass Source in dynamics
451% (matching ground truth where Source is excluded from state space)
452theta_eff(isSourceState) = 0.0;
453
454dxdt = W' * theta_eff + Alambda;
455end
456
457function theta = compute_theta(x, SQ, Sa, isSourceState)
458% COMPUTE_THETA - Compute theta vector for fluid ODE
459% For regular stations: theta = x./(SQ*x) .* min(Sa, SQ*x)
460% For Source stations (EXT scheduler): theta = 0.0
461% - Source is conceptually excluded from state space (matching ground truth)
462% - Arrivals are injected directly into queues via Alambda
463% - Setting theta = 0 prevents Source from contributing to W' * theta
464
465sumXQa = GlobalConstants.FineTol + SQ * x;
466theta = x ./ sumXQa .* min(Sa, sumXQa);
467
468% Override theta for Source stations
469% Source stations have theta = 0.0 to bypass Source in dynamics
470theta(isSourceState) = 0.0;
471end
Definition mmt.m:93