LINE Solver
MATLAB API documentation
Loading...
Searching...
No Matches
qrf_rsrd.m
1function [result, x, fval, exitflag] = qrf_rsrd(params, objective, sense)
2% QRF_RSRD - Quadratic Reduction Framework for RS-RD blocking networks
3%
4% MATLAB port of the AMPL model qrboundsrsrd_skel.mod
5%
6% Usage:
7% [result, x, fval, exitflag] = qrf_rsrd(params)
8% [result, x, fval, exitflag] = qrf_rsrd(params, objective)
9% [result, x, fval, exitflag] = qrf_rsrd(params, objective, sense)
10%
11% Inputs:
12% params - Structure with model parameters:
13% .M - Number of queues
14% .N - Total population
15% .F - [M x 1] Capacity of each queue
16% .K - [M x 1] Number of phases for each queue
17% .mu - {M x 1} cell, each mu{i} is K(i) x K(i) completion rates
18% .v - {M x 1} cell, each v{i} is K(i) x K(i) background rates
19% .r - [M x M] Routing probabilities
20% .alpha - (optional) {M x 1} cell, each alpha{i} is [N x 1] load-dependent rates
21%
22% objective - (optional) 'U1min' (default), 'U1max', or queue index 1..M
23% sense - (optional) 'min' (default) or 'max'
24%
25% Outputs:
26% result - Structure with results:
27% .U - [M x 1] Utilization of each queue
28% .Ueff - [M x 1] Effective utilization
29% .pb - [M x 1] Blocking probability
30% .p2 - Decision variable tensor (marginal probabilities)
31% x - Raw solution vector
32% fval - Objective function value
33% exitflag - Solver exit flag
34
35 if nargin < 2 || isempty(objective)
36 objective = 'U1min';
37 end
38 if nargin < 3 || isempty(sense)
39 sense = 'min';
40 end
41
42 % Extract parameters
43 M = params.M;
44 N = params.N;
45 F = params.F(:);
46 K = params.K(:);
47 mu = params.mu;
48 v = params.v;
49 r = params.r;
50 if isfield(params, 'verbose')
51 verbose = params.verbose;
52 else
53 verbose = true;
54 end
55
56 % Default alpha (load-independent)
57 if isfield(params, 'alpha') && ~isempty(params.alpha)
58 alpha = params.alpha;
59 else
60 alpha = cell(M, 1);
61 for i = 1:M
62 alpha{i} = ones(N, 1);
63 end
64 end
65
66 % Compute transition rates q(i,j,k,h,n)
67 % q{i,j} is a K(i) x K(i) x (N+1) array
68 q = cell(M, M);
69 for i = 1:M
70 for j = 1:M
71 q{i,j} = zeros(K(i), K(i), N+1);
72 for ki = 1:K(i)
73 for hi = 1:K(i)
74 for n = 1:N % n=0 gives q=0
75 if j ~= i
76 q{i,j}(ki, hi, n+1) = r(i,j) * mu{i}(ki, hi) * alpha{i}(n);
77 else
78 q{i,j}(ki, hi, n+1) = (v{i}(ki, hi) + r(i,i) * mu{i}(ki, hi)) * alpha{i}(n);
79 end
80 end
81 end
82 end
83 end
84 end
85
86 %% Build variable indexing
87 % p2(j, nj, kj, i, ni, hi) for j in 1:M, nj in 0:F(j), kj in 1:K(j),
88 % i in 1:M, ni in 0:F(i), hi in 1:K(i)
89
90 % Count variables and build index map
91 if verbose; fprintf('Building variable index map...\n'); end
92 varCount = 0;
93 p2idx = cell(M, 1);
94 for j = 1:M
95 p2idx{j} = cell(F(j)+1, K(j), M);
96 for nj = 0:F(j)
97 for kj = 1:K(j)
98 for i = 1:M
99 p2idx{j}{nj+1, kj, i} = zeros(F(i)+1, K(i));
100 for ni = 0:F(i)
101 for hi = 1:K(i)
102 varCount = varCount + 1;
103 p2idx{j}{nj+1, kj, i}(ni+1, hi) = varCount;
104 end
105 end
106 end
107 end
108 end
109 end
110 nP2Vars = varCount;
111
112 % U and Ueff variables: U(i,k,ni) and Ueff(i,k,ni) for ni >= 1
113 Uidx = cell(M, 1);
114 Ueffidx = cell(M, 1);
115 for i = 1:M
116 Uidx{i} = zeros(K(i), F(i));
117 Ueffidx{i} = zeros(K(i), F(i));
118 for ki = 1:K(i)
119 for ni = 1:F(i)
120 varCount = varCount + 1;
121 Uidx{i}(ki, ni) = varCount;
122 varCount = varCount + 1;
123 Ueffidx{i}(ki, ni) = varCount;
124 end
125 end
126 end
127
128 nVars = varCount;
129 if verbose; fprintf('Total variables: %d (p2: %d, U/Ueff: %d)\n', nVars, nP2Vars, nVars - nP2Vars); end
130
131 %% Build constraints
132 Aeq = [];
133 beq = [];
134 Aineq = [];
135 bineq = [];
136
137 % Helper to get variable index
138 getIdx = @(j, nj, kj, i, ni, hi) p2idx{j}{nj+1, kj, i}(ni+1, hi);
139 getUIdx = @(i, ki, ni) Uidx{i}(ki, ni);
140 getUeffIdx = @(i, ki, ni) Ueffidx{i}(ki, ni);
141
142 if verbose; fprintf('Building constraints...\n'); end
143
144 %% ONE: Normalization - sum over nj, kj equals 1 for each j
145 if verbose; fprintf(' ONE constraints...\n'); end
146 for j = 1:M
147 row = zeros(1, nVars);
148 for nj = 0:F(j)
149 for kj = 1:K(j)
150 idx = getIdx(j, nj, kj, j, nj, kj);
151 row(idx) = 1;
152 end
153 end
154 Aeq = [Aeq; row];
155 beq = [beq; 1];
156 end
157
158 %% ZERO constraints - fix infeasible states to zero
159 if verbose; fprintf(' ZERO constraints...\n'); end
160 lb = zeros(nVars, 1);
161 ub = ones(nVars, 1);
162 % U and Ueff bounds are [0, 1]
163 for i = 1:M
164 for ki = 1:K(i)
165 for ni = 1:F(i)
166 ub(getUIdx(i, ki, ni)) = 1;
167 ub(getUeffIdx(i, ki, ni)) = 1;
168 end
169 end
170 end
171
172 for j = 1:M
173 for nj = 0:F(j)
174 for kj = 1:K(j)
175 for i = 1:M
176 for ni = 0:F(i)
177 for hi = 1:K(i)
178 idx = getIdx(j, nj, kj, i, ni, hi);
179
180 % ZERO1: i==j, nj==ni, h<>k
181 if i == j && nj == ni && hi ~= kj
182 ub(idx) = 0;
183 end
184
185 % ZERO2: i==j, nj<>ni
186 if i == j && nj ~= ni
187 ub(idx) = 0;
188 end
189
190 % ZERO3: i<>j, nj+ni > N
191 if i ~= j && nj + ni > N
192 ub(idx) = 0;
193 end
194
195 % ZERO6: i<>j, N-nj-ni > sum of other capacities
196 if i ~= j
197 sumOtherF = sum(F) - F(i) - F(j);
198 if N - nj - ni > sumOtherF
199 ub(idx) = 0;
200 end
201 end
202 end
203 end
204 end
205 end
206 end
207
208 % ZERO7: N-nj > sum of other capacities (for diagonal)
209 for nj = 0:F(j)
210 for kj = 1:K(j)
211 sumOtherF = sum(F) - F(j);
212 if N - nj > sumOtherF
213 idx = getIdx(j, nj, kj, j, nj, kj);
214 ub(idx) = 0;
215 end
216 end
217 end
218 end
219
220 %% SYMMETRY: p2(i,ni,hi,j,nj,kj) = p2(j,nj,kj,i,ni,hi)
221 if verbose; fprintf(' SYMMETRY constraints...\n'); end
222 for j = 1:M
223 for nj = 0:F(j)
224 for kj = 1:K(j)
225 for i = 1:M
226 if i <= j
227 continue; % avoid duplicate constraints
228 end
229 for ni = 0:F(i)
230 for hi = 1:K(i)
231 idx1 = getIdx(j, nj, kj, i, ni, hi);
232 idx2 = getIdx(i, ni, hi, j, nj, kj);
233 if idx1 ~= idx2
234 row = zeros(1, nVars);
235 row(idx1) = 1;
236 row(idx2) = -1;
237 Aeq = [Aeq; row];
238 beq = [beq; 0];
239 end
240 end
241 end
242 end
243 end
244 end
245 end
246
247 %% MARGINALS: p2(j,nj,kj,j,nj,kj) = sum over ni,hi of p2(j,nj,kj,i,ni,hi)
248 if verbose; fprintf(' MARGINALS constraints...\n'); end
249 for j = 1:M
250 for kj = 1:K(j)
251 for nj = 0:F(j)
252 for i = 1:M
253 if i == j
254 continue;
255 end
256 row = zeros(1, nVars);
257 idx_diag = getIdx(j, nj, kj, j, nj, kj);
258 row(idx_diag) = 1;
259 for ni = 0:F(i)
260 for hi = 1:K(i)
261 idx = getIdx(j, nj, kj, i, ni, hi);
262 row(idx) = row(idx) - 1;
263 end
264 end
265 Aeq = [Aeq; row];
266 beq = [beq; 0];
267 end
268 end
269 end
270 end
271
272 %% UCLASSIC: U(i,k,ni) = p2(i,ni,k,i,ni,k)
273 if verbose; fprintf(' UCLASSIC constraints...\n'); end
274 for i = 1:M
275 for ki = 1:K(i)
276 for ni = 1:F(i)
277 row = zeros(1, nVars);
278 idx_U = getUIdx(i, ki, ni);
279 idx_p2 = getIdx(i, ni, ki, i, ni, ki);
280 row(idx_U) = 1;
281 row(idx_p2) = -1;
282 Aeq = [Aeq; row];
283 beq = [beq; 0];
284 end
285 end
286 end
287
288 %% UEFFS: Ueff(i,k,ni) = p2(i,ni,k,i,ni,k) - sum_{j: r(i,j)>0} r(i,j)*p2(i,ni,k,j,F(j),h)
289 if verbose; fprintf(' UEFFS constraints...\n'); end
290 for i = 1:M
291 for ki = 1:K(i)
292 for ni = 1:F(i)
293 row = zeros(1, nVars);
294 idx_Ueff = getUeffIdx(i, ki, ni);
295 idx_p2_diag = getIdx(i, ni, ki, i, ni, ki);
296 row(idx_Ueff) = 1;
297 row(idx_p2_diag) = -1;
298 % Add blocking terms
299 for j = 1:M
300 if j ~= i && r(i,j) > 0
301 for hj = 1:K(j)
302 idx_block = getIdx(i, ni, ki, j, F(j), hj);
303 row(idx_block) = r(i,j);
304 end
305 end
306 end
307 Aeq = [Aeq; row];
308 beq = [beq; 0];
309 end
310 end
311 end
312
313 %% THM2: Phase balance (Theorem 1 - THM:sdeffective)
314 % sum_{ni} (sum_{j<>i, h<>k} q(i,j,k,h,ni)*Ueff(i,k,ni) + sum_{h<>k} q(i,i,k,h,ni)*U(i,k,ni))
315 % = sum_{ni} (sum_{j<>i, h<>k} q(i,j,h,k,ni)*Ueff(i,h,ni) + sum_{h<>k} q(i,i,h,k,ni)*U(i,h,ni))
316 if verbose; fprintf(' THM2 (Phase balance) constraints...\n'); end
317 for i = 1:M
318 for ki = 1:K(i)
319 row = zeros(1, nVars);
320 % LHS terms
321 for ni = 1:F(i)
322 % Terms with Ueff (j<>i)
323 for j = 1:M
324 if j == i
325 continue;
326 end
327 for hi = 1:K(i)
328 if hi == ki
329 continue;
330 end
331 coef = q{i,j}(ki, hi, ni+1);
332 idx_Ueff = getUeffIdx(i, ki, ni);
333 row(idx_Ueff) = row(idx_Ueff) + coef;
334 end
335 end
336 % Terms with U (j==i, h<>k)
337 for hi = 1:K(i)
338 if hi == ki
339 continue;
340 end
341 coef = q{i,i}(ki, hi, ni+1);
342 idx_p2 = getIdx(i, ni, ki, i, ni, ki);
343 row(idx_p2) = row(idx_p2) + coef;
344 end
345 end
346 % RHS terms (subtract)
347 for ni = 1:F(i)
348 % Terms with Ueff (j<>i)
349 for j = 1:M
350 if j == i
351 continue;
352 end
353 for hi = 1:K(i)
354 if hi == ki
355 continue;
356 end
357 coef = q{i,j}(hi, ki, ni+1);
358 idx_Ueff = getUeffIdx(i, hi, ni);
359 row(idx_Ueff) = row(idx_Ueff) - coef;
360 end
361 end
362 % Terms with U (j==i, h<>k)
363 for hi = 1:K(i)
364 if hi == ki
365 continue;
366 end
367 coef = q{i,i}(hi, ki, ni+1);
368 idx_p2 = getIdx(i, ni, hi, i, ni, hi);
369 row(idx_p2) = row(idx_p2) - coef;
370 end
371 end
372 Aeq = [Aeq; row];
373 beq = [beq; 0];
374 end
375 end
376
377 %% THM1: Population constraint (Theorem 2)
378 % sum_i sum_ni sum_hi ni * p2(j,nj,kj,i,ni,hi) = N * p2(j,nj,kj,j,nj,kj)
379 if verbose; fprintf(' THM1 (Population) constraints...\n'); end
380 for j = 1:M
381 for kj = 1:K(j)
382 for nj = 0:F(j)
383 row = zeros(1, nVars);
384 % RHS: -N * p2(j,nj,kj,j,nj,kj)
385 idx_diag = getIdx(j, nj, kj, j, nj, kj);
386 row(idx_diag) = -N;
387 % LHS: sum
388 for i = 1:M
389 for ni = 1:F(i) % ni >= 1
390 for hi = 1:K(i)
391 idx = getIdx(j, nj, kj, i, ni, hi);
392 row(idx) = row(idx) + ni;
393 end
394 end
395 end
396 Aeq = [Aeq; row];
397 beq = [beq; 0];
398 end
399 end
400 end
401
402 %% THM3a: Marginal balance for ni in 1:F(i)-1
403 if verbose; fprintf(' THM3a (Marginal balance) constraints...\n'); end
404 for i = 1:M
405 for ni = 1:(F(i)-1)
406 row = zeros(1, nVars);
407 % LHS: arrivals to queue i
408 for j = 1:M
409 if j == i
410 continue;
411 end
412 for kj = 1:K(j)
413 for hj = 1:K(j)
414 for ui = 1:K(i)
415 for nj = 1:F(j)
416 idx = getIdx(j, nj, kj, i, ni, ui);
417 row(idx) = row(idx) + q{j,i}(kj, hj, nj+1);
418 end
419 end
420 end
421 end
422 end
423 % RHS: departures from queue i at ni+1
424 for j = 1:M
425 if j == i
426 continue;
427 end
428 for ki = 1:K(i)
429 for hi = 1:K(i)
430 for uj = 1:K(j)
431 for nj = 0:(F(j)-1)
432 idx = getIdx(i, ni+1, ki, j, nj, uj);
433 row(idx) = row(idx) - q{i,j}(ki, hi, ni+2);
434 end
435 end
436 end
437 end
438 end
439 Aeq = [Aeq; row];
440 beq = [beq; 0];
441 end
442 end
443
444 %% THM3b: Marginal balance for ni=0, per phase
445 if verbose; fprintf(' THM3b (Marginal balance ni=0) constraints...\n'); end
446 for i = 1:M
447 for ui = 1:K(i)
448 row = zeros(1, nVars);
449 % LHS: arrivals to queue i at ni=0
450 for j = 1:M
451 if j == i
452 continue;
453 end
454 for kj = 1:K(j)
455 for hj = 1:K(j)
456 for nj = 1:F(j)
457 idx = getIdx(j, nj, kj, i, 0, ui);
458 row(idx) = row(idx) + q{j,i}(kj, hj, nj+1);
459 end
460 end
461 end
462 end
463 % RHS: departures from queue i at ni=1
464 for j = 1:M
465 if j == i
466 continue;
467 end
468 for ki = 1:K(i)
469 for nj = 0:(F(j)-1)
470 for hj = 1:K(j)
471 idx = getIdx(i, 1, ki, j, nj, hj);
472 row(idx) = row(idx) - q{i,j}(ki, ui, 2);
473 end
474 end
475 end
476 end
477 Aeq = [Aeq; row];
478 beq = [beq; 0];
479 end
480 end
481
482 %% QBAL: Queue balance constraint
483 % This is a complex balance equation that tightens the bounds
484 if verbose; fprintf(' QBAL (Queue balance) constraints...\n'); end
485 for i = 1:M
486 for ki = 1:K(i)
487 row = zeros(1, nVars);
488
489 % LHS Term 1: sum{h<>k} sum{j<>i} sum{ni} sum{u} sum{nj} q[i,j,k,h,ni]*ni*p2[i,ni,k,j,nj,u]
490 for hi = 1:K(i)
491 if hi == ki
492 continue;
493 end
494 for j = 1:M
495 if j == i
496 continue;
497 end
498 for ni = 1:F(i)
499 for uj = 1:K(j)
500 for nj = 0:(F(j)-1)
501 coef = q{i,j}(ki, hi, ni+1) * ni;
502 idx = getIdx(i, ni, ki, j, nj, uj);
503 row(idx) = row(idx) + coef;
504 end
505 end
506 end
507 end
508 end
509
510 % LHS Term 2: sum{h<>k} sum{ni} q[i,i,k,h,ni]*ni*p2[i,ni,k,i,ni,k]
511 for hi = 1:K(i)
512 if hi == ki
513 continue;
514 end
515 for ni = 1:F(i)
516 coef = q{i,i}(ki, hi, ni+1) * ni;
517 idx = getIdx(i, ni, ki, i, ni, ki);
518 row(idx) = row(idx) + coef;
519 end
520 end
521
522 % LHS Term 3: sum{j<>i} sum{h} sum{ni} sum{u} sum{nj<=min(F(j)-1,N-ni)} q[i,j,h,k,ni]*p2[i,ni,h,j,nj,u]
523 for j = 1:M
524 if j == i
525 continue;
526 end
527 for hi = 1:K(i)
528 for ni = 1:F(i)
529 for uj = 1:K(j)
530 maxNj = min(F(j)-1, N-ni);
531 for nj = 0:maxNj
532 coef = q{i,j}(hi, ki, ni+1);
533 idx = getIdx(i, ni, hi, j, nj, uj);
534 row(idx) = row(idx) + coef;
535 end
536 end
537 end
538 end
539 end
540
541 % RHS Term 1: sum{j<>i} sum{h} sum{ni<=F(i)-1} sum{u} sum{nj>=1} q[j,i,h,u,nj]*p2[i,ni,k,j,nj,h]
542 for j = 1:M
543 if j == i
544 continue;
545 end
546 for hj = 1:K(j)
547 for ni = 0:(F(i)-1)
548 for uj = 1:K(j)
549 for nj = 1:F(j)
550 coef = q{j,i}(hj, uj, nj+1);
551 idx = getIdx(i, ni, ki, j, nj, hj);
552 row(idx) = row(idx) - coef;
553 end
554 end
555 end
556 end
557 end
558
559 % RHS Term 2: sum{h<>k} sum{ni} q[i,i,h,k,ni]*ni*p2[i,ni,h,i,ni,h]
560 for hi = 1:K(i)
561 if hi == ki
562 continue;
563 end
564 for ni = 1:F(i)
565 coef = q{i,i}(hi, ki, ni+1) * ni;
566 idx = getIdx(i, ni, hi, i, ni, hi);
567 row(idx) = row(idx) - coef;
568 end
569 end
570
571 % RHS Term 3: sum{h<>k} sum{j<>i} sum{ni} sum{u} sum{nj} q[i,j,h,k,ni]*ni*p2[i,ni,h,j,nj,u]
572 for hi = 1:K(i)
573 if hi == ki
574 continue;
575 end
576 for j = 1:M
577 if j == i
578 continue;
579 end
580 for ni = 1:F(i)
581 for uj = 1:K(j)
582 for nj = 0:(F(j)-1)
583 coef = q{i,j}(hi, ki, ni+1) * ni;
584 idx = getIdx(i, ni, hi, j, nj, uj);
585 row(idx) = row(idx) - coef;
586 end
587 end
588 end
589 end
590 end
591
592 Aeq = [Aeq; row];
593 beq = [beq; 0];
594 end
595 end
596
597 %% THM4: Queue-length bound inequality (Theorem 5)
598 % sum_t sum_ht sum_nj sum_nt nt*p2(j,nj,kj,t,nt,ht) >= N * sum_hi sum_nj sum_ni p2(j,nj,kj,i,ni,hi)
599 if verbose; fprintf(' THM4 (Queue-length bound) constraints...\n'); end
600 for j = 1:M
601 for kj = 1:K(j)
602 for i = 1:M
603 row = zeros(1, nVars);
604 % LHS: sum_t sum_ht sum_nj sum_nt nt * p2
605 for t = 1:M
606 for ht = 1:K(t)
607 for nj = 0:F(j)
608 for nt = 1:F(t) % nt >= 1
609 idx = getIdx(j, nj, kj, t, nt, ht);
610 row(idx) = row(idx) + nt;
611 end
612 end
613 end
614 end
615 % RHS: -N * sum
616 for hi = 1:K(i)
617 for nj = 0:F(j)
618 for ni = 1:F(i) % ni >= 1
619 idx = getIdx(j, nj, kj, i, ni, hi);
620 row(idx) = row(idx) - N;
621 end
622 end
623 end
624 Aineq = [Aineq; -row]; % >= becomes <= with negation
625 bineq = [bineq; 0];
626 end
627 end
628 end
629
630 %% Build objective function
631 if verbose; fprintf('Building objective function...\n'); end
632 c = zeros(nVars, 1);
633
634 if ischar(objective)
635 if strcmp(objective, 'U1min') || strcmp(objective, 'U1max')
636 targetQueue = 1;
637 else
638 error('Unknown objective: %s', objective);
639 end
640 else
641 targetQueue = objective;
642 end
643
644 % Utilization = sum over k, n of p2(i,n,k,i,n,k)
645 for ki = 1:K(targetQueue)
646 for ni = 1:F(targetQueue)
647 idx = getIdx(targetQueue, ni, ki, targetQueue, ni, ki);
648 c(idx) = 1;
649 end
650 end
651
652 if strcmp(sense, 'max')
653 c = -c;
654 end
655
656 %% Solve LP
657 if verbose; fprintf('Solving LP with %d variables and %d equality + %d inequality constraints...\n', ...
658 nVars, size(Aeq, 1), size(Aineq, 1)); end
659
660 if verbose
661 options = optimoptions('linprog', 'Display', 'final', 'Algorithm', 'interior-point');
662 else
663 options = optimoptions('linprog', 'Display', 'off', 'Algorithm', 'interior-point');
664 end
665
666 [x, fval, exitflag] = linprog(c, Aineq, bineq, Aeq, beq, lb, ub, options);
667
668 if strcmp(sense, 'max')
669 fval = -fval;
670 end
671
672 %% Extract results
673 result = struct();
674 result.objective = fval;
675 result.exitflag = exitflag;
676
677 if exitflag > 0
678 % Compute utilizations from U and Ueff variables
679 result.U = zeros(M, 1);
680 result.Ueff = zeros(M, 1);
681 result.pb = zeros(M, 1);
682
683 for i = 1:M
684 % U(i) = sum over k, ni of U(i,k,ni)
685 for ki = 1:K(i)
686 for ni = 1:F(i)
687 idx_U = getUIdx(i, ki, ni);
688 idx_Ueff = getUeffIdx(i, ki, ni);
689 result.U(i) = result.U(i) + x(idx_U);
690 result.Ueff(i) = result.Ueff(i) + x(idx_Ueff);
691 end
692 end
693 result.pb(i) = result.U(i) - result.Ueff(i);
694 end
695
696 % Store p2 as a function handle for easy access
697 result.getP2 = @(j, nj, kj, i, ni, hi) x(getIdx(j, nj, kj, i, ni, hi));
698 end
699
700 if verbose; fprintf('\n=== Results ===\n'); end
701 if verbose; fprintf('Objective value: %f\n', fval); end
702 if verbose; fprintf('Exit flag: %d\n', exitflag); end
703 if exitflag > 0 && verbose
704 fprintf('\nUtilizations:\n');
705 for i = 1:M
706 fprintf(' Queue %d: U = %.6f, Ueff = %.6f, pb = %.6f\n', ...
707 i, result.U(i), result.Ueff(i), result.pb(i));
708 end
709 end
710end