def MMAPPH1NPPR(D, sigma, S, *argv): """ Returns various performane measures of a continuous time MMAP[K]/PH[K]/1 non-preemptive priority queue, see [1]_. Parameters ---------- D : list of matrices of shape (N,N), length (K+1) The D0...DK matrices of the arrival process. D1 corresponds to the lowest, DK to the highest priority. sigma : list of row vectors, length (K) The list containing the initial probability vectors of the service time distributions of the various customer types. The length of the vectors does not have to be the same. S : list of square matrices, length (K) The transient generators of the phase type distributions representing the service time of the jobs belonging to various types. further parameters : The rest of the function parameters specify the options and the performance measures to be computed. The supported performance measures and options in this function are: +----------------+--------------------+----------------------------------------+ | Parameter name | Input parameters | Output | +================+====================+========================================+ | "ncMoms" | Number of moments | The moments of the number of customers | +----------------+--------------------+----------------------------------------+ | "ncDistr" | Upper limit K | The distribution of the number of | | | | customers from level 0 to level K-1 | +----------------+--------------------+----------------------------------------+ | "stMoms" | Number of moments | The sojourn time moments | +----------------+--------------------+----------------------------------------+ | "stDistr" | A vector of points | The sojourn time distribution at the | | | | requested points (cummulative, cdf) | +----------------+--------------------+----------------------------------------+ | "prec" | The precision | Numerical precision used as a stopping | | | | condition when solving the Riccati and | | | | the matrix-quadratic equations | +----------------+--------------------+----------------------------------------+ | "erlMaxOrder" | Integer number | The maximal Erlang order used in the | | | | erlangization procedure. The default | | | | value is 200. | +----------------+--------------------+----------------------------------------+ | "classes" | Vector of integers | Only the performance measures | | | | belonging to these classes are | | | | returned. If not given, all classes | | | | are analyzed. | +----------------+--------------------+----------------------------------------+ (The quantities related to the number of customers in the system include the customer in the server, and the sojourn time related quantities include the service times as well) Returns ------- Ret : list of the performance measures Each entry of the list corresponds to a performance measure requested. Each entry is a matrix, where the columns belong to the various job types. If there is just a single item, then it is not put into a list. References ---------- .. [1] G. Horvath, "Efficient analysis of the MMAP[K]/PH[K]/1 priority queue", European Journal of Operational Research, 246(1), 128-139, 2015. """ K = len(D) - 1 # parse options eaten = [] erlMaxOrder = 200 precision = 1e-14 classes = np.arange(0, K) for i in range(len(argv)): if argv[i] == "prec": precision = argv[i + 1] eaten.append(i) eaten.append(i + 1) elif argv[i] == "erlMaxOrder": erlMaxOrder = argv[i + 1] eaten.append(i) eaten.append(i + 1) elif argv[i] == "classes": classes = np.array(argv[i + 1]) - 1 eaten.append(i) eaten.append(i + 1) if butools.checkInput and not CheckMMAPRepresentation(D): raise Exception( 'MMAPPH1PRPR: The arrival process is not a valid MMAP representation!' ) if butools.checkInput: for k in range(K): if not CheckPHRepresentation(sigma[k], S[k]): raise Exception( 'MMAPPH1PRPR: the vector and matrix describing the service times is not a valid PH representation!' ) # some preparation D0 = D[0] N = D0.shape[0] I = ml.eye(N) sD = ml.zeros((N, N)) for Di in D: sD += Di s = [] M = np.empty(K) for i in range(K): s.append(np.sum(-S[i], 1)) M[i] = sigma[i].size # step 1. solution of the workload process of the joint queue # =========================================================== sM = np.sum(M) Qwmm = ml.matrix(D0) Qwpm = ml.zeros((N * sM, N)) Qwmp = ml.zeros((N, N * sM)) Qwpp = ml.zeros((N * sM, N * sM)) kix = 0 for i in range(K): Qwmp[:, kix:kix + N * M[i]] = np.kron(D[i + 1], sigma[i]) Qwpm[kix:kix + N * M[i], :] = np.kron(I, s[i]) Qwpp[kix:kix + N * M[i], :][:, kix:kix + N * M[i]] = np.kron(I, S[i]) kix += N * M[i] # calculate fundamental matrices Psiw, Kw, Uw = FluidFundamentalMatrices(Qwpp, Qwpm, Qwmp, Qwmm, 'PKU', precision) # calculate boundary vector Ua = ml.ones((N, 1)) + 2 * np.sum(Qwmp * (-Kw).I, 1) pm = Linsolve( ml.hstack((Uw, Ua)).T, ml.hstack((ml.zeros((1, N)), ml.ones((1, 1)))).T).T ro = ((1.0 - np.sum(pm)) / 2.0) / ( np.sum(pm) + (1.0 - np.sum(pm)) / 2.0 ) # calc idle time with weight=1, and the busy time with weight=1/2 kappa = pm / np.sum(pm) pi = CTMCSolve(sD) lambd = [] for i in range(K): lambd.append(np.sum(pi * D[i + 1])) Psiw = [] Qwmp = [] Qwzp = [] Qwpp = [] Qwmz = [] Qwpz = [] Qwzz = [] Qwmm = [] Qwpm = [] Qwzm = [] for k in range(K): # step 2. construct a workload process for classes k...K # ====================================================== Mlo = np.sum(M[:k]) Mhi = np.sum(M[k:]) Qkwpp = ml.zeros((N * Mlo * Mhi + N * Mhi, N * Mlo * Mhi + N * Mhi)) Qkwpz = ml.zeros((N * Mlo * Mhi + N * Mhi, N * Mlo)) Qkwpm = ml.zeros((N * Mlo * Mhi + N * Mhi, N)) Qkwmz = ml.zeros((N, N * Mlo)) Qkwmp = ml.zeros((N, N * Mlo * Mhi + N * Mhi)) Dlo = ml.matrix(D0) for i in range(k): Dlo = Dlo + D[i + 1] Qkwmm = Dlo Qkwzp = ml.zeros((N * Mlo, N * Mlo * Mhi + N * Mhi)) Qkwzm = ml.zeros((N * Mlo, N)) Qkwzz = ml.zeros((N * Mlo, N * Mlo)) kix = 0 for i in range(k, K): kix2 = 0 for j in range(k): bs = N * M[j] * M[i] bs2 = N * M[j] Qkwpp[kix:kix + bs, kix:kix + bs] = np.kron(I, np.kron(ml.eye(M[j]), S[i])) Qkwpz[kix:kix + bs, kix2:kix2 + bs2] = np.kron(I, np.kron(ml.eye(M[j]), s[i])) Qkwzp[kix2:kix2 + bs2, kix:kix + bs] = np.kron(D[i + 1], np.kron(ml.eye(M[j]), sigma[i])) kix += bs kix2 += bs2 for i in range(k, K): bs = N * M[i] Qkwpp[kix:kix + bs, :][:, kix:kix + bs] = np.kron(I, S[i]) Qkwpm[kix:kix + bs, :] = np.kron(I, s[i]) Qkwmp[:, kix:kix + bs] = np.kron(D[i + 1], sigma[i]) kix += bs kix = 0 for j in range(k): bs = N * M[j] Qkwzz[kix:kix + bs, kix:kix + bs] = np.kron(Dlo, ml.eye(M[j])) + np.kron(I, S[j]) Qkwzm[kix:kix + bs, :] = np.kron(I, s[j]) kix += bs if Qkwzz.shape[0] > 0: Psikw = FluidFundamentalMatrices( Qkwpp + Qkwpz * (-Qkwzz).I * Qkwzp, Qkwpm + Qkwpz * (-Qkwzz).I * Qkwzm, Qkwmp, Qkwmm, 'P', precision) else: Psikw = FluidFundamentalMatrices(Qkwpp, Qkwpm, Qkwmp, Qkwmm, 'P', precision) Psiw.append(Psikw) Qwzp.append(Qkwzp) Qwmp.append(Qkwmp) Qwpp.append(Qkwpp) Qwmz.append(Qkwmz) Qwpz.append(Qkwpz) Qwzz.append(Qkwzz) Qwmm.append(Qkwmm) Qwpm.append(Qkwpm) Qwzm.append(Qkwzm) # step 3. calculate Phi vectors # ============================= lambdaS = sum(lambd) phi = [(1 - ro) * kappa * (-D0) / lambdaS] q0 = [[]] qL = [[]] for k in range(K - 1): sDk = ml.matrix(D0) for j in range(k + 1): sDk = sDk + D[j + 1] # pk pk = sum(lambd[:k + 1]) / lambdaS - (1 - ro) * kappa * np.sum( sDk, 1) / lambdaS # A^(k,1) Qwzpk = Qwzp[k + 1] vix = 0 Ak = [] for ii in range(k + 1): bs = N * M[ii] V1 = Qwzpk[vix:vix + bs, :] Ak.append( np.kron(I, sigma[ii]) * (-np.kron(sDk, ml.eye(M[ii])) - np.kron(I, S[ii])).I * (np.kron(I, s[ii]) + V1 * Psiw[k + 1])) vix += bs # B^k Qwmpk = Qwmp[k + 1] Bk = Qwmpk * Psiw[k + 1] ztag = phi[0] * ((-D0).I * D[k + 1] * Ak[k] - Ak[0] + (-D0).I * Bk) for i in range(k): ztag += phi[i + 1] * (Ak[i] - Ak[i + 1]) + phi[0] * ( -D0).I * D[i + 1] * Ak[i] Mx = ml.eye(Ak[k].shape[0]) - Ak[k] Mx[:, 0] = ml.ones((N, 1)) phi.append( ml.hstack((pk, ztag[:, 1:])) * Mx.I) # phi(k) = Psi^(k)_k * p(k). Psi^(k)_i = phi(i) / p(k) q0.append(phi[0] * (-D0).I) qLii = [] for ii in range(k + 1): qLii.append((phi[ii + 1] - phi[ii] + phi[0] * (-D0).I * D[ii + 1]) * np.kron(I, sigma[ii]) * (-np.kron(sDk, ml.eye(M[ii])) - np.kron(I, S[ii])).I) qL.append(ml.hstack(qLii)) # step 4. calculate performance measures # ====================================== Ret = [] for k in classes: sD0k = ml.matrix(D0) for i in range(k): sD0k += D[i + 1] if k < K - 1: # step 4.1 calculate distribution of the workload process right # before the arrivals of class k jobs # ============================================================ if Qwzz[k].shape[0] > 0: Kw = Qwpp[k] + Qwpz[k] * ( -Qwzz[k]).I * Qwzp[k] + Psiw[k] * Qwmp[k] else: Kw = Qwpp[k] + Psiw[k] * Qwmp[k] BM = ml.zeros((0, 0)) CM = ml.zeros((0, N)) DM = ml.zeros((0, 0)) for i in range(k): BM = la.block_diag(BM, np.kron(I, S[i])) CM = ml.vstack((CM, np.kron(I, s[i]))) DM = la.block_diag(DM, np.kron(D[k + 1], ml.eye(M[i]))) if k > 0: Kwu = ml.vstack((ml.hstack( (Kw, (Qwpz[k] + Psiw[k] * Qwmz[k]) * (-Qwzz[k]).I * DM)), ml.hstack((ml.zeros( (BM.shape[0], Kw.shape[1])), BM)))) Bwu = ml.vstack((Psiw[k] * D[k + 1], CM)) iniw = ml.hstack( (q0[k] * Qwmp[k] + qL[k] * Qwzp[k], qL[k] * DM)) pwu = q0[k] * D[k + 1] else: Kwu = Kw Bwu = Psiw[k] * D[k + 1] iniw = pm * Qwmp[k] pwu = pm * D[k + 1] norm = np.sum(pwu) + np.sum(iniw * (-Kwu).I * Bwu) pwu = pwu / norm iniw = iniw / norm # step 4.2 create the fluid model whose first passage time equals the # WAITING time of the low prioroity customers # ================================================================== KN = Kwu.shape[0] Qspp = ml.zeros( (KN + N * np.sum(M[k + 1:]), KN + N * np.sum(M[k + 1:]))) Qspm = ml.zeros((KN + N * np.sum(M[k + 1:]), N)) Qsmp = ml.zeros((N, KN + N * np.sum(M[k + 1:]))) Qsmm = sD0k + D[k + 1] kix = 0 for i in range(k + 1, K): bs = N * M[i] Qspp[KN + kix:KN + kix + bs, :][:, KN + kix:KN + kix + bs] = np.kron(I, S[i]) Qspm[KN + kix:KN + kix + bs, :] = np.kron(I, s[i]) Qsmp[:, KN + kix:KN + kix + bs] = np.kron(D[i + 1], sigma[i]) kix += bs Qspp[:KN, :][:, :KN] = Kwu Qspm[:KN, :] = Bwu inis = ml.hstack((iniw, ml.zeros((1, N * np.sum(M[k + 1:]))))) # calculate fundamental matrix Psis = FluidFundamentalMatrices(Qspp, Qspm, Qsmp, Qsmm, 'P', precision) # step 4.3. calculate the performance measures # ========================================== argIx = 0 while argIx < len(argv): if argIx in eaten: argIx += 1 continue elif type(argv[argIx]) is str and argv[argIx] == "stMoms": # MOMENTS OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfSTMoms = argv[argIx + 1] # calculate waiting time moments Pn = [Psis] wtMoms = [] for n in range(1, numOfSTMoms + 1): A = Qspp + Psis * Qsmp B = Qsmm + Qsmp * Psis C = -2 * n * Pn[n - 1] bino = 1 for i in range(1, n): bino = bino * (n - i + 1) / i C += bino * Pn[i] * Qsmp * Pn[n - i] P = la.solve_sylvester(A, B, -C) Pn.append(P) wtMoms.append(np.sum(inis * P * (-1)**n) / 2**n) # calculate RESPONSE time moments Pnr = [np.sum(inis * Pn[0]) * sigma[k]] rtMoms = [] for n in range(1, numOfSTMoms + 1): P = n * Pnr[n - 1] * (-S[k]).I + (-1)**n * np.sum( inis * Pn[n]) * sigma[k] / 2**n Pnr.append(P) rtMoms.append( np.sum(P) + np.sum(pwu) * math.factorial(n) * np.sum(sigma[k] * (-S[k]).I**n)) Ret.append(rtMoms) argIx += 1 elif type(argv[argIx]) is str and argv[argIx] == "stDistr": # DISTRIBUTION OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ stCdfPoints = argv[argIx + 1] res = [] for t in stCdfPoints: L = erlMaxOrder lambdae = L / t / 2 Psie = FluidFundamentalMatrices( Qspp - lambdae * ml.eye(Qspp.shape[0]), Qspm, Qsmp, Qsmm - lambdae * ml.eye(Qsmm.shape[0]), 'P', precision) Pn = [Psie] pr = (np.sum(pwu) + np.sum(inis * Psie)) * (1 - np.sum( sigma[k] * (ml.eye(S[k].shape[0]) - S[k] / 2 / lambdae).I**L)) for n in range(1, L): A = Qspp + Psie * Qsmp - lambdae * ml.eye( Qspp.shape[0]) B = Qsmm + Qsmp * Psie - lambdae * ml.eye( Qsmm.shape[0]) C = 2 * lambdae * Pn[n - 1] for i in range(1, n): C += Pn[i] * Qsmp * Pn[n - i] P = la.solve_sylvester(A, B, -C) Pn.append(P) pr += np.sum(inis * P) * ( 1 - np.sum(sigma[k] * (np.eye(S[k].shape[0]) - S[k] / 2 / lambdae).I**(L - n))) res.append(pr) Ret.append(np.array(res)) argIx += 1 elif type(argv[argIx]) is str and (argv[argIx] == "ncMoms" or argv[argIx] == "ncDistr"): W = (-np.kron(sD - D[k + 1], ml.eye(M[k])) - np.kron(I, S[k])).I * np.kron(D[k + 1], ml.eye(M[k])) iW = (ml.eye(W.shape[0]) - W).I w = np.kron(ml.eye(N), sigma[k]) omega = (-np.kron(sD - D[k + 1], ml.eye(M[k])) - np.kron(I, S[k])).I * np.kron(I, s[k]) if argv[argIx] == "ncMoms": # MOMENTS OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLMoms = argv[argIx + 1] # first calculate it at departure instants Psii = [Psis] QLDPn = [inis * Psii[0] * w * iW] for n in range(1, numOfQLMoms + 1): A = Qspp + Psis * Qsmp B = Qsmm + Qsmp * Psis C = n * Psii[n - 1] * D[k + 1] bino = 1 for i in range(1, n): bino = bino * (n - i + 1) / i C = C + bino * Psii[i] * Qsmp * Psii[n - i] P = la.solve_sylvester(A, B, -C) Psii.append(P) QLDPn.append(n * QLDPn[n - 1] * iW * W + inis * P * w * iW) for n in range(numOfQLMoms + 1): QLDPn[n] = (QLDPn[n] + pwu * w * iW**(n + 1) * W**n) * omega # now calculate it at random time instance QLPn = [pi] qlMoms = [] iTerm = (ml.ones((N, 1)) * pi - sD).I for n in range(1, numOfQLMoms + 1): sumP = np.sum(QLDPn[n]) + n * np.sum( (QLDPn[n - 1] - QLPn[n - 1] * D[k + 1] / lambd[k]) * iTerm * D[k + 1]) P = sumP * pi + n * ( QLPn[n - 1] * D[k + 1] - QLDPn[n - 1] * lambd[k]) * iTerm QLPn.append(P) qlMoms.append(np.sum(P)) qlMoms = MomsFromFactorialMoms(qlMoms) Ret.append(qlMoms) argIx += 1 elif argv[argIx] == "ncDistr": # DISTRIBUTION OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLProbs = argv[argIx + 1] Psid = FluidFundamentalMatrices( Qspp, Qspm, Qsmp, sD0k, 'P', precision) Pn = [Psid] XDn = inis * Psid * w dqlProbs = (XDn + pwu * w) * omega for n in range(1, numOfQLProbs): A = Qspp + Psid * Qsmp B = sD0k + Qsmp * Psid C = Pn[n - 1] * D[k + 1] for i in range(1, n): C += Pn[i] * Qsmp * Pn[n - i] P = la.solve_sylvester(A, B, -C) Pn.append(P) XDn = XDn * W + inis * P * w dqlProbs = ml.vstack( (dqlProbs, (XDn + pwu * w * W**n) * omega)) # now calculate it at random time instance iTerm = -(sD - D[k + 1]).I qlProbs = lambd[k] * dqlProbs[0, :] * iTerm for n in range(1, numOfQLProbs): P = (qlProbs[n - 1, :] * D[k + 1] + lambd[k] * (dqlProbs[n, :] - dqlProbs[n - 1, :])) * iTerm qlProbs = ml.vstack((qlProbs, P)) qlProbs = np.sum(qlProbs, 1).A.flatten() Ret.append(qlProbs) argIx += 1 else: raise Exception("MMAPPH1NPPR: Unknown parameter " + str(argv[argIx])) argIx += 1 elif k == K - 1: # step 3. calculate the performance measures # ========================================== argIx = 0 while argIx < len(argv): if argIx in eaten: argIx += 1 continue elif type(argv[argIx]) is str and (argv[argIx] == "stMoms" or argv[argIx] == "stDistr"): Kw = Qwpp[k] + Qwpz[k] * ( -Qwzz[k]).I * Qwzp[k] + Psiw[k] * Qwmp[k] AM = ml.zeros((0, 0)) BM = ml.zeros((0, 0)) CM = ml.zeros((0, 1)) DM = ml.zeros((0, 0)) for i in range(k): AM = la.block_diag( AM, np.kron(ml.ones((N, 1)), np.kron(ml.eye(M[i]), s[k]))) BM = la.block_diag(BM, S[i]) CM = ml.vstack((CM, s[i])) DM = la.block_diag(DM, np.kron(D[k + 1], ml.eye(M[i]))) Z = ml.vstack((ml.hstack( (Kw, ml.vstack((AM, ml.zeros( (N * M[k], AM.shape[1])))))), ml.hstack((ml.zeros( (BM.shape[0], Kw.shape[1])), BM)))) z = ml.vstack((ml.zeros( (AM.shape[0], 1)), np.kron(ml.ones((N, 1)), s[k]), CM)) iniw = ml.hstack((q0[k] * Qwmp[k] + qL[k] * Qwzp[k], ml.zeros((1, BM.shape[0])))) zeta = iniw / np.sum(iniw * (-Z).I * z) if argv[argIx] == "stMoms": # MOMENTS OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfSTMoms = argv[argIx + 1] rtMoms = [] for i in range(1, numOfSTMoms + 1): rtMoms.append( np.sum( math.factorial(i) * zeta * (-Z).I**(i + 1) * z)) Ret.append(rtMoms) argIx += 1 if argv[argIx] == "stDistr": # DISTRIBUTION OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ stCdfPoints = argv[argIx + 1] rtDistr = [] for t in stCdfPoints: rtDistr.append( np.sum(zeta * (-Z).I * (ml.eye(Z.shape[0]) - la.expm(Z * t)) * z)) Ret.append(np.array(rtDistr)) argIx += 1 elif type(argv[argIx]) is str and (argv[argIx] == "ncMoms" or argv[argIx] == "ncDistr"): L = ml.zeros((N * np.sum(M), N * np.sum(M))) B = ml.zeros((N * np.sum(M), N * np.sum(M))) F = ml.zeros((N * np.sum(M), N * np.sum(M))) kix = 0 for i in range(K): bs = N * M[i] F[kix:kix + bs, :][:, kix:kix + bs] = np.kron( D[k + 1], ml.eye(M[i])) L[kix:kix + bs, :][:, kix:kix + bs] = np.kron( sD0k, ml.eye(M[i])) + np.kron(I, S[i]) if i < K - 1: L[kix:kix + bs, :][:, N * np.sum(M[:k]):] = np.kron( I, s[i] * sigma[k]) else: B[kix:kix + bs, :][:, N * np.sum(M[:k]):] = np.kron( I, s[i] * sigma[k]) kix += bs R = QBDFundamentalMatrices(B, L, F, 'R', precision) p0 = ml.hstack((qL[k], q0[k] * np.kron(I, sigma[k]))) p0 = p0 / np.sum(p0 * (ml.eye(R.shape[0]) - R).I) if argv[argIx] == "ncMoms": # MOMENTS OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLMoms = argv[argIx + 1] qlMoms = [] for i in range(1, numOfQLMoms + 1): qlMoms.append( np.sum( math.factorial(i) * p0 * R**i * (ml.eye(R.shape[0]) - R).I**(i + 1))) Ret.append(MomsFromFactorialMoms(qlMoms)) elif argv[argIx] == "ncDistr": # DISTRIBUTION OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLProbs = argv[argIx + 1] qlProbs = [np.sum(p0)] for i in range(1, numOfQLProbs): qlProbs.append(np.sum(p0 * R**i)) Ret.append(np.array(qlProbs)) argIx += 1 else: raise Exception("MMAPPH1NPPR: Unknown parameter " + str(argv[argIx])) argIx += 1 if len(Ret) == 1: return Ret[0] else: return Ret
def MAPMAP1(D0, D1, S0, S1, *argv): """ Returns various performane measures of a continuous time MAP/MAP/1 queue. In a MAP/MAP/1 queue both the arrival and the service processes are characterized by Markovian arrival processes. Parameters ---------- D0 : matrix, shape(N,N) The transitions of the arrival MAP not accompanied by job arrivals D1 : matrix, shape(N,N) The transitions of the arrival MAP accompanied by job arrivals S0 : matrix, shape(N,N) The transitions of the service MAP not accompanied by job service S1 : matrix, shape(N,N) The transitions of the service MAP accompanied by job service further parameters : The rest of the function parameters specify the options and the performance measures to be computed. The supported performance measures and options in this function are: +----------------+--------------------+----------------------------------------+ | Parameter name | Input parameters | Output | +================+====================+========================================+ | "ncMoms" | Number of moments | The moments of the number of customers | +----------------+--------------------+----------------------------------------+ | "ncDistr" | Upper limit K | The distribution of the number of | | | | customers from level 0 to level K-1 | +----------------+--------------------+----------------------------------------+ | "ncDistrMG" | None | The vector-matrix parameters of the | | | | matrix-geometric distribution of the | | | | number of customers in the system | +----------------+--------------------+----------------------------------------+ | "ncDistrDPH" | None | The vector-matrix parameters of the | | | | matrix-geometric distribution of the | | | | number of customers in the system, | | | | converted to a discrete PH | | | | representation | +----------------+--------------------+----------------------------------------+ | "stMoms" | Number of moments | The sojourn time moments | +----------------+--------------------+----------------------------------------+ | "stDistr" | A vector of points | The sojourn time distribution at the | | | | requested points (cummulative, cdf) | +----------------+--------------------+----------------------------------------+ | "stDistrME" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | sojourn time distribution | +----------------+--------------------+----------------------------------------+ | "stDistrPH" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | sojourn time distribution, converted | | | | to a continuous PH representation | +----------------+--------------------+----------------------------------------+ | "prec" | The precision | Numerical precision used as a stopping | | | | condition when solving the | | | | matrix-quadratic equation | +----------------+--------------------+----------------------------------------+ (The quantities related to the number of customers in the system include the customer in the server, and the sojourn time related quantities include the service times as well) Returns ------- Ret : list of the performance measures Each entry of the list corresponds to a performance measure requested. If there is just a single item, then it is not put into a list. Notes ----- "ncDistrMG" and "stDistrME" behave much better numerically than "ncDistrDPH" and "stDistrPH". """ # parse options prec = 1e-14 needST = False eaten = [] for i in range(len(argv)): if argv[i] == "prec": prec = argv[i + 1] eaten.append(i) eaten.append(i + 1) elif type(argv[i]) is str and len( argv[i]) > 2 and argv[i][0:2] == "st": needST = True if butools.checkInput and not CheckMAPRepresentation(D0, D1): raise Exception( 'MAPMAP1: The arrival process (D0,D1) is not a valid MAP representation!' ) if butools.checkInput and not CheckMAPRepresentation(S0, S1): raise Exception( 'MAPMAP1: The service process (S0,S1) is not a valid MAP representation!' ) IA = ml.eye(D0.shape[0]) IS = ml.eye(S0.shape[0]) B = np.kron(IA, S1) L = np.kron(D0, IS) + np.kron(IA, S0) F = np.kron(D1, IS) L0 = np.kron(D0, IS) pi0, R = QBDSolve(B, L, F, L0, prec) N = pi0.shape[1] I = ml.eye(N) if needST: # calculate the distribution of the age at departures U = L + R * B Rh = -U.I * F T = np.kron(IA, S0) + Rh * B eta = pi0 * F * (I - Rh).I eta = eta / np.sum(eta) Ret = [] argIx = 0 while argIx < len(argv): if argIx in eaten: argIx += 1 continue elif type(argv[argIx]) is str and argv[argIx] == "ncDistrDPH": # transform it to DPH alpha = pi0 * R * (I - R).I A = Diag(alpha).I * R.T * Diag(alpha) Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx] == "ncDistrMG": # transform it to MG B = SimilarityMatrixForVectors(np.sum((I - R).I * R, 1), np.ones((N, 1))) Bi = B.I A = B * R * Bi alpha = pi0 * Bi Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx] == "ncMoms": numOfMoms = argv[argIx + 1] argIx += 1 moms = [] iR = (I - R).I for m in range(1, numOfMoms + 1): moms.append( math.factorial(m) * np.sum(pi0 * iR**(m + 1) * R**m)) Ret.append(MomsFromFactorialMoms(moms)) elif type(argv[argIx]) is str and argv[argIx] == "ncDistr": numOfQLProbs = argv[argIx + 1] argIx += 1 values = np.empty(numOfQLProbs) values[0] = np.sum(pi0) RPow = I for p in range(numOfQLProbs - 1): RPow = RPow * R values[p + 1] = np.sum(pi0 * RPow) Ret.append(values) elif type(argv[argIx]) is str and argv[argIx] == "stDistrPH": # transform it to PH representation beta = CTMCSolve(S0 + S1) theta = DTMCSolve(-D0.I * D1) vv = np.kron(theta, beta) ix = np.arange(N) nz = ix[vv.flat > prec] delta = Diag(vv[:, nz]) alpha = ml.ones((1, N)) * B[nz, :].T * delta / np.sum(beta * S1) A = delta.I * T[nz, :][:, nz].T * delta Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx] == "stDistrME": Ret.append(eta) Ret.append(T) elif type(argv[argIx]) is str and argv[argIx] == "stMoms": numOfMoms = argv[argIx + 1] argIx += 1 moms = [] iT = -T.I for m in range(1, numOfMoms + 1): moms.append(math.factorial(m) * np.sum(eta * iT**m)) Ret.append(moms) elif type(argv[argIx]) is str and argv[argIx] == "stDistr": points = argv[argIx + 1] argIx += 1 values = np.empty(points.shape) for p in range(len(points.flat)): values.flat[p] = 1.0 - np.sum( eta * la.expm(T * points.flat[p])) Ret.append(values) else: raise Exception("MAPMAP1: Unknown parameter " + str(argv[argIx])) argIx += 1 if len(Ret) == 1: return Ret[0] else: return Ret
def RandomPH(order, mean=1.0, zeroEntries=0, maxTrials=1000, prec=1e-7): """ Returns a random phase-type distribution with a given order. Parameters ---------- order : int The size of the phase-type distribution mean : double, optional The mean of the phase-type distribution zeroEntries : int, optional The number of zero entries in the initial vector, generator matrix and closing vector maxTrials : int, optional The maximum number of trials to find a proper PH (that has an irreducible phase process and none of its parameters is all-zero). The default value is 1000. prec : double, optional Numerical precision for checking the irreducibility. The default value is 1e-14. Returns ------- alpha : vector, shape (1,M) The initial probability vector of the phase-type distribution. A : matrix, shape (M,M) The transient generator matrix of the phase-type distribution. Notes ----- If the procedure fails, try to increase the 'maxTrials' parameter. """ if zeroEntries > (order + 1) * (order - 1): raise Exception( "RandomPH: Too many zero entries requested! Try to decrease the zero entries number!" ) # distribute the zero entries among the rows def allZeroDistr(states, zeros): if states == 1: return [[zeros]] else: o = [] for i in range(zeros + 1): x = allZeroDistr(states - 1, zeros - i) for j in range(len(x)): xt = x[j] xt.append(i) xt.sort() # check if we have it already if o.count(xt) == 0: o.append(xt) return o zeroDistr = allZeroDistr(order, zeroEntries) trials = 1 while trials < maxTrials: # select a configuration from zeroDistr: it is a list describing the zero entries in each row zdix = np.random.permutation(len(zeroDistr)) for k in range(len(zeroDistr)): zDistr = zeroDistr[zdix[k]] B = np.zeros((order, order + 2)) for i in range(order): rp = np.random.permutation(order + 1) a = np.zeros(order + 1) for j in range(order + 1 - zDistr[i]): a[rp[j]] = np.random.rand() B[i, 0:i] = a[0:i] B[i, i + 1:] = a[i:] # construct PH parameters A = ml.matrix(B[:, :order]) a = ml.matrix(B[:, order + 1]).T A = A - Diag(np.sum(A, 1) + a) alpha = ml.matrix(B[:, order]) # check if it is a proper PH (irreducible phase process & no full zero matrix) if np.all(A == 0.0) or np.all(alpha == 0.0) or np.all(a == 0.0): continue alpha = alpha / np.sum(alpha) D = A + a * alpha if la.matrix_rank(D) == order - 1: pi = CTMCSolve(D) if np.min(np.abs(pi)) > prec: # scale to the mean value m = MomentsFromPH(alpha, A, 1)[0] A *= m / mean return (alpha, A) trials += 1 raise Exception("No feasible random PH found with such many zero entries!")
def FluFluQueue(Qin, Rin, Qout, Rout, srv0stop, *argv): """ Returns various performane measures of a fluid queue with independent fluid arrival and service processes. Two types of boundary behavior is available. If srv0stop=false, the output process evolves continuously even if the queue is empty. If srv0stop=true, the output process slows down if there is fewer fluid in the queue than it can serve. If the queue is empty and the fluid input rate is zero, the output process freezes till fluid arrives. Parameters ---------- Qin : matrix, shape (N,N) The generator of the background Markov chain corresponding to the input process Rin : matrix, shape (N,N) Diagonal matrix containing the fluid input rates associated to the states of the input background process Qout : matrix, shape (N,N) The generator of the background Markov chain corresponding to the output process Rout : matrix, shape (N,N) Diagonal matrix containing the fluid output rates associated to the states of the input background process srv0stop : bool If true, the service output process slows down if there is fewer fluid in the queue than it can serve. If false, the output process evolves continuously. further parameters : The rest of the function parameters specify the options and the performance measures to be computed. The supported performance measures and options in this function are: +----------------+--------------------+--------------------------------------+ | Parameter name | Input parameters | Output | +================+====================+======================================+ | "flMoms" | Number of moments | The moments of the fluid level | +----------------+--------------------+--------------------------------------+ | "flDistr" | A vector of points | The fluid level distribution at | | | | the requested points (cdf) | +----------------+--------------------+--------------------------------------+ | "flDistrME" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | fluid level distribution | +----------------+--------------------+--------------------------------------+ | "flDistrPH" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | fluid level distribution, converted | | | | to a PH representation | +----------------+--------------------+--------------------------------------+ | "stMoms" | Number of moments | The sojourn time moments of fluid | | | | drops | +----------------+--------------------+--------------------------------------+ | "stDistr" | A vector of points | The sojourn time distribution at the | | | | requested points (cummulative, cdf) | +----------------+--------------------+--------------------------------------+ | "stDistrME" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | sojourn time distribution | +----------------+--------------------+--------------------------------------+ | "stDistrPH" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | sojourn time distribution, converted | | | | to a PH representation | +----------------+--------------------+--------------------------------------+ | "prec" | The precision | Numerical precision to check if the | | | | input is valid and it is also used | | | | as a stopping condition when solving | | | | the Riccati equation | +----------------+--------------------+--------------------------------------+ Returns ------- Ret : list of the performance measures Each entry of the list corresponds to a performance measure requested. If there is just a single item, then it is not put into a list. Notes ----- "flDistrME" and "stDistrME" behave much better numerically than "flDistrPH" and "stDistrPH". References ---------- .. [1] Horvath G, Telek M, "Sojourn times in fluid queues with independent and dependent input and output processes PERFORMANCE EVALUATION 79: pp. 160-181, 2014. """ # parse options prec = 1e-14 needST = False needQL = False Q0 = [] eaten = [] for i in range(len(argv)): if type(argv[i]) is str and argv[i]=="prec": prec = argv[i+1] eaten.append(i) eaten.append(i+1) elif type(argv[i]) is str and len(argv[i])>2 and argv[i][0:2]=="st": needST = True elif type(argv[i]) is str and len(argv[i])>2 and argv[i][0:2]=="fl": needQL = True if butools.checkInput and not CheckGenerator(Qin,False): raise Exception('FluFluQueue: Generator matrix Qin is not Markovian!') if butools.checkInput and not CheckGenerator(Qout,False): raise Exception('FluFluQueue: Generator matrix Qout is not Markovian!') if butools.checkInput and (np.any(np.diag(Rin)<-butools.checkPrecision) or np.any(np.diag(Rout)<-butools.checkPrecision)): raise Exception('FluFluQueue: Fluid rates Rin and Rout must be non-negative !') Iin = ml.eye(Qin.shape[0]) Iout = ml.eye(Qout.shape[0]) if needQL: Q = np.kron(Qin,Iout)+np.kron(Iin,Qout) if srv0stop: Q0 = np.kron(Qin,Iout)+np.kron(Rin, la.pinv(Rout)*Qout) else: Q0 = Q mass0, ini, K, clo = GeneralFluidSolve (Q, np.kron(Rin,Iout)-np.kron(Iin,Rout), Q0, prec) if needST: Rh = np.kron(Rin,Iout) - np.kron(Iin,Rout) Qh = np.kron(Qin, Rout) + np.kron(Rin, Qout) massh, inih, Kh, cloh = GeneralFluidSolve (Qh, Rh, prec=prec) # sojourn time density in case of # srv0stop = false: inih*expm(Kh*x)*cloh*kron(Rin,Iout)/lambda # srv0stop = true: inih*expm(Kh*x)*cloh*kron(Rin,Rout)/lambda/mu lambd = np.sum(CTMCSolve(Qin)*Rin) mu = np.sum(CTMCSolve(Qout)*Rout) Ret = [] argIx = 0 while argIx<len(argv): if argIx in eaten: argIx += 1 continue elif type(argv[argIx]) is str and argv[argIx]=='flDistrPH': # transform it to PH Delta = Diag(Linsolve(K.T,-ini.T)) # Delta = diag (ini*inv(-K)); A = Delta.I*K.T*Delta alpha = np.sum(clo,1).T*Delta Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx]=='flDistrME': # transform it to ME B = SimilarityMatrixForVectors((-K).I*np.sum(clo,1), np.ones((K.shape[0],1))) Bi = B.I alpha = ini*Bi A = B*K*Bi Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx]=='flMoms': numOfMoms = argv[argIx+1] argIx += 1 moms = [] iK = -K.I for m in range(1,numOfMoms+1): moms.append(math.factorial(m)*np.sum(ini*iK**(m+1)*clo)) Ret.append(moms) elif type(argv[argIx]) is str and argv[argIx]=='flDistr': points = argv[argIx+1] argIx += 1 values = np.empty(points.shape) iK = -K.I for p in range(len(points.flat)): values.flat[p] = np.sum(mass0) + np.sum(ini*(ml.eye(K.shape[0])-la.expm(K*points.flat[p]))*iK*clo) Ret.append (values) elif type(argv[argIx]) is str and argv[argIx]=='stDistrPH': # convert result to PH representation Delta = Diag(Linsolve(Kh.T,-inih.T)) # Delta = diag (inih*inv(-Kh)); A = Delta.I*Kh.T*Delta if not srv0stop: alpha = np.sum(Delta*cloh*np.kron(Rin,Iout)/lambd,1).T else: alpha = np.sum(Delta*cloh*np.kron(Rin,Rout)/lambd/mu,1).T Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx]=='stDistrME': # convert result to ME representation if not srv0stop: B = SimilarityMatrixForVectors(np.sum(cloh*np.kron(Rin,Iout)/lambd,1), np.ones((Kh.shape[0],1))) else: B = SimilarityMatrixForVectors(np.sum(cloh*np.kron(Rin,Rout)/lambd/mu,1), np.ones((Kh.shape[0],1))) iB = B.I A = B*Kh*iB alpha = inih*(-Kh).I*iB Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx]=='stMoms': numOfMoms = argv[argIx+1] argIx += 1 moms = [] if srv0stop: kclo = cloh*np.kron(Rin,Rout)/lambd/mu else: kclo = cloh*np.kron(Rin,Iout)/lambd iKh = -Kh.I for m in range(1,numOfMoms+1): moms.append(math.factorial(m)*np.sum(inih*iKh**(m+1)*kclo)) Ret.append(moms) elif type(argv[argIx]) is str and argv[argIx]=='stDistr': points = argv[argIx+1] argIx += 1 values = np.empty(points.shape) if srv0stop: kclo = cloh*np.kron(Rin,Rout)/lambd/mu else: kclo = cloh*np.kron(Rin,Iout)/lambd iKh = -Kh.I for p in range(len(points.flat)): values.flat[p] = 1.0-np.sum(inih*la.expm(Kh*points.flat[p])*iKh*kclo) Ret.append(values) else: raise Exception ("FluFluQueue: Unknown parameter "+str(argv[argIx])) argIx += 1 if len(Ret)==1: return Ret[0] else: return Ret
def MMAPPH1FCFS(D, sigma, S, *argv): """ Returns various performane measures of a MMAP[K]/PH[K]/1 first-come-first-serve queue, see [1]_. Parameters ---------- D : list of matrices of shape (N,N), length (K+1) The D0...DK matrices of the arrival process. sigma : list of row vectors, length (K) The list containing the initial probability vectors of the service time distributions of the various customer types. The length of the vectors does not have to be the same. S : list of square matrices, length (K) The transient generators of the phase type distributions representing the service time of the jobs belonging to various types. further parameters : The rest of the function parameters specify the options and the performance measures to be computed. The supported performance measures and options in this function are: +----------------+--------------------+----------------------------------------+ | Parameter name | Input parameters | Output | +================+====================+========================================+ | "ncMoms" | Number of moments | The moments of the number of customers | +----------------+--------------------+----------------------------------------+ | "ncDistr" | Upper limit K | The distribution of the number of | | | | customers from level 0 to level K-1 | +----------------+--------------------+----------------------------------------+ | "stMoms" | Number of moments | The sojourn time moments | +----------------+--------------------+----------------------------------------+ | "stDistr" | A vector of points | The sojourn time distribution at the | | | | requested points (cummulative, cdf) | +----------------+--------------------+----------------------------------------+ | "stDistrME" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | sojourn time distribution | +----------------+--------------------+----------------------------------------+ | "stDistrPH" | None | The vector-matrix parameters of the | | | | matrix-exponentially distributed | | | | sojourn time distribution, converted | | | | to a continuous PH representation | +----------------+--------------------+----------------------------------------+ | "prec" | The precision | Numerical precision used as a stopping | | | | condition when solving the Riccati | | | | equation | +----------------+--------------------+----------------------------------------+ | "classes" | Vector of integers | Only the performance measures | | | | belonging to these classes are | | | | returned. If not given, all classes | | | | are analyzed. | +----------------+--------------------+----------------------------------------+ (The quantities related to the number of customers in the system include the customer in the server, and the sojourn time related quantities include the service times as well) Returns ------- Ret : list of the performance measures Each entry of the list corresponds to a performance measure requested. Each entry is a matrix, where the columns belong to the various job types. If there is just a single item, then it is not put into a list. References ---------- .. [1] Qiming He, "Analysis of a continuous time SM[K]/PH[K]/1/FCFS queue: Age process, sojourn times, and queue lengths", Journal of Systems Science and Complexity, 25(1), pp 133-155, 2012. """ K = len(D) - 1 # parse options eaten = [] precision = 1e-14 classes = np.arange(0, K) for i in range(len(argv)): if argv[i] == "prec": precision = argv[i + 1] eaten.append(i) eaten.append(i + 1) elif argv[i] == "classes": classes = np.array(argv[i + 1]) - 1 eaten.append(i) eaten.append(i + 1) if butools.checkInput and not CheckMMAPRepresentation(D): raise Exception( 'MMAPPH1FCFS: The arrival process is not a valid MMAP representation!' ) if butools.checkInput: for k in range(K): if not CheckPHRepresentation(sigma[k], S[k]): raise Exception( 'MMAPPH1FCFS: the vector and matrix describing the service times is not a valid PH representation!' ) # some preparation D0 = D[0] N = D0.shape[0] Ia = ml.eye(N) Da = ml.zeros((N, N)) for q in range(K): Da += D[q + 1] theta = CTMCSolve(D0 + Da) beta = [CTMCSolve(S[k] + ml.sum(-S[k], 1) * sigma[k]) for k in range(K)] lambd = [np.sum(theta * D[k + 1]) for k in range(K)] mu = [np.sum(beta[k] * (-S[k])) for k in range(K)] Nsk = [S[k].shape[0] for k in range(K)] ro = np.sum(np.array(lambd) / np.array(mu)) alpha = theta * Da / sum(lambd) D0i = (-D0).I Sa = S[0] sa = [ml.zeros(sigma[0].shape)] * K sa[0] = sigma[0] ba = [ml.zeros(beta[0].shape)] * K ba[0] = beta[0] sv = [ml.zeros((Nsk[0], 1))] * K sv[0] = ml.sum(-S[0], 1) Pk = [D0i * D[q + 1] for q in range(K)] for k in range(1, K): Sa = la.block_diag(Sa, S[k]) for q in range(K): if q == k: sa[q] = ml.hstack((sa[q], sigma[k])) ba[q] = ml.hstack((ba[q], beta[k])) sv[q] = ml.vstack((sv[q], -np.sum(S[k], 1))) else: sa[q] = ml.hstack((sa[q], ml.zeros(sigma[k].shape))) ba[q] = ml.hstack((ba[q], ml.zeros(beta[k].shape))) sv[q] = ml.vstack((sv[q], ml.zeros((Nsk[k], 1)))) Sa = ml.matrix(Sa) P = D0i * Da iVec = ml.kron(D[1], sa[0]) for k in range(1, K): iVec += ml.kron(D[k + 1], sa[k]) Ns = Sa.shape[0] Is = ml.eye(Ns) # step 1. solve the age process of the queue # ========================================== # solve Y0 and calculate T Y0 = FluidFundamentalMatrices(ml.kron(Ia, Sa), ml.kron(Ia, -ml.sum(Sa, 1)), iVec, D0, "P", precision) T = ml.kron(Ia, Sa) + Y0 * iVec # calculate pi0 and v0 pi0 = ml.zeros((1, T.shape[0])) for k in range(K): pi0 += ml.kron(theta * D[k + 1], ba[k] / mu[k]) pi0 = -pi0 * T iT = (-T).I oa = ml.ones((N, 1)) # step 2. calculate performance measures # ====================================== Ret = [] for k in classes: argIx = 0 clo = iT * ml.kron(oa, sv[k]) while argIx < len(argv): if argIx in eaten: argIx += 1 continue elif type(argv[argIx]) is str and argv[argIx] == "stMoms": numOfSTMoms = argv[argIx + 1] rtMoms = [] for m in range(1, numOfSTMoms + 1): rtMoms.append( math.factorial(m) * np.sum(pi0 * iT**m * clo / (pi0 * clo))) Ret.append(rtMoms) argIx += 1 elif type(argv[argIx]) is str and argv[argIx] == "stDistr": stCdfPoints = argv[argIx + 1] cdf = [] for t in stCdfPoints: pr = 1 - np.sum(pi0 * la.expm(T * t) * clo / (pi0 * clo)) cdf.append(pr) Ret.append(np.array(cdf)) argIx += 1 elif type(argv[argIx]) is str and argv[argIx] == "stDistrME": Bm = SimilarityMatrixForVectors(clo / (pi0 * clo), ml.ones((N * Ns, 1))) Bmi = Bm.I A = Bm * T * Bmi alpha = pi0 * Bmi Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx] == "stDistrPH": vv = pi0 * iT ix = np.arange(N * Ns) nz = ix[vv.flat > precision] delta = Diag(vv[:, nz]) cl = -T * clo / (pi0 * clo) alpha = cl[nz, :].T * delta A = delta.I * T[nz, :][:, nz].T * delta Ret.append(alpha) Ret.append(A) elif type(argv[argIx]) is str and argv[argIx] == "ncDistr": numOfQLProbs = argv[argIx + 1] argIx += 1 values = np.empty(numOfQLProbs) jm = ml.zeros((Ns, 1)) jm[np.sum(Nsk[0:k]):np.sum(Nsk[0:k + 1]), :] = 1 jmc = ml.ones((Ns, 1)) jmc[np.sum(Nsk[0:k]):np.sum(Nsk[0:k + 1]), :] = 0 LmCurr = la.solve_sylvester(T, ml.kron(D0 + Da - D[k + 1], Is), -ml.eye(N * Ns)) values[0] = 1 - ro + np.sum(pi0 * LmCurr * ml.kron(oa, jmc)) for i in range(1, numOfQLProbs): LmPrev = LmCurr LmCurr = la.solve_sylvester( T, ml.kron(D0 + Da - D[k + 1], Is), -LmPrev * ml.kron(D[k + 1], Is)) values[i] = np.sum(pi0 * LmCurr * ml.kron(oa, jmc) + pi0 * LmPrev * ml.kron(oa, jm)) Ret.append(values) elif type(argv[argIx]) is str and argv[argIx] == "ncMoms": numOfQLMoms = argv[argIx + 1] argIx += 1 jm = ml.zeros((Ns, 1)) jm[np.sum(Nsk[0:k]):np.sum(Nsk[0:k + 1]), :] = 1 ELn = [ la.solve_sylvester(T, ml.kron(D0 + Da, Is), -ml.eye(N * Ns)) ] qlMoms = [] for n in range(1, numOfQLMoms + 1): bino = 1 Btag = ml.zeros((N * Ns, N * Ns)) for i in range(n): Btag += bino * ELn[i] bino *= (n - i) / (i + 1) ELn.append( la.solve_sylvester(T, ml.kron(D0 + Da, Is), -Btag * ml.kron(D[k + 1], Is))) qlMoms.append( np.sum(pi0 * ELn[n]) + np.sum(pi0 * Btag * ml.kron(oa, jm))) Ret.append(qlMoms) else: raise Exception("MMAPPH1FCFS: Unknown parameter " + str(argv[argIx])) argIx += 1 if len(Ret) == 1: return Ret[0] else: return Ret
def MMAPPH1PRPR(D, sigma, S, *argv): """ Returns various performane measures of a MMAP[K]/PH[K]/1 preemptive resume priority queue, see [1]_. Parameters ---------- D : list of matrices of shape (N,N), length (K+1) The D0...DK matrices of the arrival process. D1 corresponds to the lowest, DK to the highest priority. sigma : list of row vectors, length (K) The list containing the initial probability vectors of the service time distributions of the various customer types. The length of the vectors does not have to be the same. S : list of square matrices, length (K) The transient generators of the phase type distributions representing the service time of the jobs belonging to various types. further parameters : The rest of the function parameters specify the options and the performance measures to be computed. The supported performance measures and options in this function are: +----------------+--------------------+----------------------------------------+ | Parameter name | Input parameters | Output | +================+====================+========================================+ | "ncMoms" | Number of moments | The moments of the number of customers | +----------------+--------------------+----------------------------------------+ | "ncDistr" | Upper limit K | The distribution of the number of | | | | customers from level 0 to level K-1 | +----------------+--------------------+----------------------------------------+ | "stMoms" | Number of moments | The sojourn time moments | +----------------+--------------------+----------------------------------------+ | "stDistr" | A vector of points | The sojourn time distribution at the | | | | requested points (cummulative, cdf) | +----------------+--------------------+----------------------------------------+ | "prec" | The precision | Numerical precision used as a stopping | | | | condition when solving the Riccati and | | | | the matrix-quadratic equations | +----------------+--------------------+----------------------------------------+ | "erlMaxOrder" | Integer number | The maximal Erlang order used in the | | | | erlangization procedure. The default | | | | value is 200. | +----------------+--------------------+----------------------------------------+ | "classes" | Vector of integers | Only the performance measures | | | | belonging to these classes are | | | | returned. If not given, all classes | | | | are analyzed. | +----------------+--------------------+----------------------------------------+ (The quantities related to the number of customers in the system include the customer in the server, and the sojourn time related quantities include the service times as well) Returns ------- Ret : list of the performance measures Each entry of the list corresponds to a performance measure requested. Each entry is a matrix, where the columns belong to the various job types. If there is just a single item, then it is not put into a list. References ---------- .. [1] G. Horvath, "Efficient analysis of the MMAP[K]/PH[K]/1 priority queue", European Journal of Operational Research, 246(1), 128-139, 2015. """ K = len(D) - 1 # parse options eaten = [] erlMaxOrder = 200 precision = 1e-14 classes = np.arange(0, K) for i in range(len(argv)): if argv[i] == "prec": precision = argv[i + 1] eaten.append(i) eaten.append(i + 1) elif argv[i] == "erlMaxOrder": erlMaxOrder = argv[i + 1] eaten.append(i) eaten.append(i + 1) elif argv[i] == "classes": classes = np.array(argv[i + 1]) - 1 eaten.append(i) eaten.append(i + 1) if butools.checkInput and not CheckMMAPRepresentation(D): raise Exception( 'MMAPPH1PRPR: The arrival process is not a valid MMAP representation!' ) if butools.checkInput: for k in range(K): if not CheckPHRepresentation(sigma[k], S[k]): raise Exception( 'MMAPPH1PRPR: the vector and matrix describing the service times is not a valid PH representation!' ) # some preparation D0 = D[0] N = D0.shape[0] I = ml.eye(N) sD = ml.zeros((N, N)) for Di in D: sD += Di s = [] M = np.empty(K) for i in range(K): s.append(np.sum(-S[i], 1)) M[i] = sigma[i].size Ret = [] for k in classes: # step 1. solution of the workload process of the system # ====================================================== sM = np.sum(M[k:K]) Qwmm = ml.matrix(D0) for i in range(k): Qwmm += D[i + 1] Qwpm = ml.zeros((N * sM, N)) Qwmp = ml.zeros((N, N * sM)) Qwpp = ml.zeros((N * sM, N * sM)) kix = 0 for i in range(k, K): Qwmp[:, kix:kix + N * M[i]] = np.kron(D[i + 1], sigma[i]) Qwpm[kix:kix + N * M[i], :] = np.kron(I, s[i]) Qwpp[kix:kix + N * M[i], :][:, kix:kix + N * M[i]] = np.kron(I, S[i]) kix += N * M[i] # calculate fundamental matrices Psiw, Kw, Uw = FluidFundamentalMatrices(Qwpp, Qwpm, Qwmp, Qwmm, 'PKU', precision) # calculate boundary vector Ua = ml.ones((N, 1)) + 2 * np.sum(Qwmp * (-Kw).I, 1) pm = Linsolve( ml.hstack((Uw, Ua)).T, ml.hstack((ml.zeros((1, N)), ml.ones((1, 1)))).T).T Bw = ml.zeros((N * sM, N)) Bw[0:N * M[k], :] = np.kron(I, s[k]) kappa = pm * Qwmp / np.sum(pm * Qwmp * (-Kw).I * Bw) if k < K - 1: # step 2. construct fluid model for the remaining sojourn time process # ==================================================================== # (for each class except the highest priority) Qsmm = ml.matrix(D0) for i in range(k + 1): Qsmm += D[i + 1] Np = Kw.shape[0] Qspm = ml.zeros((Np + N * np.sum(M[k + 1:]), N)) Qsmp = ml.zeros((N, Np + N * np.sum(M[k + 1:]))) Qspp = ml.zeros( (Np + N * np.sum(M[k + 1:]), Np + N * np.sum(M[k + 1:]))) Qspp[:Np, :Np] = Kw Qspm[:Np, :N] = Bw kix = Np for i in range(k + 1, K): Qsmp[:, kix:kix + N * M[i]] = np.kron(D[i + 1], sigma[i]) Qspm[kix:kix + N * M[i], :] = np.kron(I, s[i]) Qspp[kix:kix + N * M[i], kix:kix + N * M[i]] = np.kron(I, S[i]) kix += N * M[i] inis = ml.hstack((kappa, ml.zeros((1, N * np.sum(M[k + 1:]))))) Psis = FluidFundamentalMatrices(Qspp, Qspm, Qsmp, Qsmm, 'P', precision) # step 3. calculate the performance measures # ========================================== argIx = 0 while argIx < len(argv): if argIx in eaten: argIx += 1 continue elif type(argv[argIx]) is str and argv[argIx] == "stMoms": # MOMENTS OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfSTMoms = argv[argIx + 1] Pn = [Psis] rtMoms = [] for n in range(1, numOfSTMoms + 1): A = Qspp + Psis * Qsmp B = Qsmm + Qsmp * Psis C = -2 * n * Pn[n - 1] bino = 1 for i in range(1, n): bino = bino * (n - i + 1) / i C += bino * Pn[i] * Qsmp * Pn[n - i] P = la.solve_sylvester(A, B, -C) Pn.append(P) rtMoms.append(np.sum(inis * P * (-1)**n) / 2**n) Ret.append(rtMoms) argIx += 1 elif type(argv[argIx]) is str and argv[argIx] == "stDistr": # DISTRIBUTION OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ stCdfPoints = argv[argIx + 1] res = [] for t in stCdfPoints: L = erlMaxOrder lambd = L / t / 2 Psie = FluidFundamentalMatrices( Qspp - lambd * ml.eye(Qspp.shape[0]), Qspm, Qsmp, Qsmm - lambd * ml.eye(Qsmm.shape[0]), 'P', precision) Pn = [Psie] pr = np.sum(inis * Psie) for n in range(1, L): A = Qspp + Psie * Qsmp - lambd * ml.eye( Qspp.shape[0]) B = Qsmm + Qsmp * Psie - lambd * ml.eye( Qsmm.shape[0]) C = 2 * lambd * Pn[n - 1] for i in range(1, n): C += Pn[i] * Qsmp * Pn[n - i] P = la.solve_sylvester(A, B, -C) Pn.append(P) pr += np.sum(inis * P) res.append(pr) Ret.append(np.array(res)) argIx += 1 elif type(argv[argIx]) is str and argv[argIx] == "ncMoms": # MOMENTS OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLMoms = argv[argIx + 1] # first calculate it at departure instants QLDPn = [Psis] dqlMoms = [] for n in range(1, numOfQLMoms + 1): A = Qspp + Psis * Qsmp B = Qsmm + Qsmp * Psis C = n * QLDPn[n - 1] * D[k + 1] bino = 1 for i in range(1, n): bino = bino * (n - i + 1) / i C = C + bino * QLDPn[i] * Qsmp * QLDPn[n - i] P = la.solve_sylvester(A, B, -C) QLDPn.append(P) dqlMoms.append(np.sum(inis * P)) dqlMoms = MomsFromFactorialMoms(dqlMoms) # now calculate it at random time instance pi = CTMCSolve(sD) lambdak = np.sum(pi * D[k + 1]) QLPn = [pi] qlMoms = [] iTerm = (ml.ones((N, 1)) * pi - sD).I for n in range(1, numOfQLMoms + 1): sumP = np.sum(inis * QLDPn[n]) + n * ( inis * QLDPn[n - 1] - QLPn[n - 1] * D[k + 1] / lambdak) * iTerm * np.sum(D[k + 1], 1) P = sumP * pi + n * (QLPn[n - 1] * D[k + 1] - inis * QLDPn[n - 1] * lambdak) * iTerm QLPn.append(P) qlMoms.append(np.sum(P)) qlMoms = MomsFromFactorialMoms(qlMoms) Ret.append(qlMoms) argIx += 1 elif type(argv[argIx]) is str and argv[argIx] == "ncDistr": # DISTRIBUTION OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLProbs = argv[argIx + 1] sDk = ml.matrix(D0) for i in range(k): sDk += D[i + 1] # first calculate it at departure instants Psid = FluidFundamentalMatrices(Qspp, Qspm, Qsmp, sDk, 'P', precision) Pn = [Psid] dqlProbs = inis * Psid for n in range(1, numOfQLProbs): A = Qspp + Psid * Qsmp B = sDk + Qsmp * Psid C = Pn[n - 1] * D[k + 1] for i in range(1, n): C += Pn[i] * Qsmp * Pn[n - i] P = la.solve_sylvester(A, B, -C) Pn.append(P) dqlProbs = ml.vstack((dqlProbs, inis * P)) # now calculate it at random time instance pi = CTMCSolve(sD) lambdak = np.sum(pi * D[k + 1]) iTerm = -(sD - D[k + 1]).I qlProbs = lambdak * dqlProbs[0, :] * iTerm for n in range(1, numOfQLProbs): P = (qlProbs[n - 1, :] * D[k + 1] + lambdak * (dqlProbs[n, :] - dqlProbs[n - 1, :])) * iTerm qlProbs = ml.vstack((qlProbs, P)) qlProbs = np.sum(qlProbs, 1).A.flatten() Ret.append(qlProbs) argIx += 1 else: raise Exception("MMAPPH1PRPR: Unknown parameter " + str(argv[argIx])) argIx += 1 elif k == K - 1: # step 3. calculate the performance measures # ========================================== argIx = 0 while argIx < len(argv): if argIx in eaten: argIx += 1 continue elif type(argv[argIx]) is str and argv[argIx] == "stMoms": # MOMENTS OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfSTMoms = argv[argIx + 1] rtMoms = [] for i in range(1, numOfSTMoms + 1): rtMoms.append( np.sum( math.factorial(i) * kappa * (-Kw).I**(i + 1) * Bw)) Ret.append(rtMoms) argIx += 1 elif type(argv[argIx]) is str and argv[argIx] == "stDistr": # DISTRIBUTION OF THE SOJOURN TIME # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ stCdfPoints = argv[argIx + 1] rtDistr = [] for t in stCdfPoints: rtDistr.append( np.sum(kappa * (-Kw).I * (ml.eye(Kw.shape[0]) - la.expm(Kw * t)) * Bw)) Ret.append(np.array(rtDistr)) argIx += 1 elif type(argv[argIx]) is str and (argv[argIx] == "ncMoms" or argv[argIx] == "ncDistr"): L = np.kron(sD - D[k + 1], ml.eye(M[k])) + np.kron( ml.eye(N), S[k]) B = np.kron(ml.eye(N), s[k] * sigma[k]) F = np.kron(D[k + 1], ml.eye(M[k])) L0 = np.kron(sD - D[k + 1], ml.eye(M[k])) R = QBDFundamentalMatrices(B, L, F, 'R', precision) p0 = CTMCSolve(L0 + R * B) p0 = p0 / np.sum(p0 * (ml.eye(R.shape[0]) - R).I) if argv[argIx] == "ncMoms": # MOMENTS OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLMoms = argv[argIx + 1] qlMoms = [] for i in range(1, numOfQLMoms + 1): qlMoms.append( np.sum( math.factorial(i) * p0 * R**i * (ml.eye(R.shape[0]) - R).I**(i + 1))) Ret.append(MomsFromFactorialMoms(qlMoms)) elif argv[argIx] == "ncDistr": # DISTRIBUTION OF THE NUMBER OF JOBS # ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ numOfQLProbs = argv[argIx + 1] qlProbs = [np.sum(p0)] for i in range(1, numOfQLProbs): qlProbs.append(np.sum(p0 * R**i)) Ret.append(np.array(qlProbs)) argIx += 1 else: raise Exception("MMAPPH1PRPR: Unknown parameter " + str(argv[argIx])) argIx += 1 if len(Ret) == 1: return Ret[0] else: return Ret
def SamplesFromMMAP(D, k, initial=None, prec=1e-14): """ Generates random samples from a marked Markovian arrival process. Parameters ---------- D : list of matrices of shape(M,M), length(N) The D0...DN matrices of the MMAP K : integer The number of samples to generate. prec : double, optional Numerical precision to check if the input MMAP is valid. The default value is 1e-14. Returns ------- x : matrix, shape(K,2) The random samples. Each row consists of two columns: the inter-arrival time and the type of the arrival. """ if butools.checkInput and not CheckMMAPRepresentation(D): raise Exception( "SamplesFromMMAP: Input is not a valid MMAP representation!") N = D[0].shape[0] if initial == None: # draw initial state according to the stationary distribution stst = CTMCSolve(SumMatrixList(D)).A.flatten() cummInitial = np.cumsum(stst) r = rand() state = 0 while cummInitial[state] <= r: state += 1 else: state = initial # auxilary variables sojourn = -1.0 / np.diag(D[0]) nextpr = ml.matrix(np.diag(sojourn)) * D[0] nextpr = nextpr - ml.matrix(np.diag(np.diag(nextpr))) for i in range(1, len(D)): nextpr = np.hstack((nextpr, np.diag(sojourn) * D[i])) nextpr = np.cumsum(nextpr, 1) if len(D) > 2: x = np.empty((k, 2)) else: x = np.empty(k) for n in range(k): time = 0 # play state transitions while state < N: time -= np.log(rand()) * sojourn[state] r = rand() nstate = 0 while nextpr[state, nstate] <= r: nstate += 1 state = nstate if len(D) > 2: x[n, 0] = time x[n, 1] = state // N else: x[n] = time state = state % N return x
def RandomMMAP(order, types, mean=1.0, zeroEntries=0, maxTrials=1000, prec=1e-7): """ Returns a random Markovian arrival process with given mean value. Parameters ---------- order : int The size of the MAP types : int The number of different arrival types mean : double, optional The mean inter-arrival times of the MMAP zeroEntries : int, optional The number of zero entries in the D0 and D1 matrices maxTrials : int, optional The maximum number of trials to find a proper MMAP (that has an irreducible phase process and none of its parameters is all-zero) prec : double, optional Numerical precision for checking the irreducibility. The default value is 1e-14. Returns ------- D : list/cell of matrices of shape(M,M), length(types+1) The D0...Dtypes matrices of the MMAP """ # distribute the zero entries among the rows def allZeroDistr(states, zeros): if states == 1: return [[zeros]] else: o = [] for i in range(zeros + 1): x = allZeroDistr(states - 1, zeros - i) for j in range(len(x)): xt = x[j] xt.append(i) xt.sort() # check if we have it already if o.count(xt) == 0: o.append(xt) return o if zeroEntries > (types + 1) * order * order - 2 * order: raise Exception( "RandomMAP/MMAP: Too many zero entries requested! Try to decrease the zeroEntries parameter!" ) zeroDistr = allZeroDistr(order, zeroEntries) trials = 1 while trials < maxTrials: # select a configuration from zeroDistr: it is a list describing the zero entries in each row zdix = np.random.permutation(len(zeroDistr)) for k in range(len(zeroDistr)): zDistr = zeroDistr[zdix[k]] bad = False for d in zDistr: if d >= (types + 1) * order - 1: bad = True break if bad: continue B = np.zeros((order, (types + 1) * order)) for i in range(order): rp = np.random.permutation((types + 1) * order - 1) a = np.zeros((types + 1) * order - 1) for j in range((types + 1) * order - 1 - zDistr[i]): a[rp[j]] = np.random.rand() B[i, 0:i] = a[0:i] B[i, i + 1:] = a[i:] # construct MMAP matrices D = [] for i in range(types + 1): Di = ml.matrix(B[:, i * order:(i + 1) * order]) D.append(Di) D[0] -= Diag(np.sum(B, 1)) # check if it is a proper MAP (irreducible phase process & no full zero matrix) sumD = ml.zeros((order, order)) for i in range(len(D)): sumD += D[i] if la.matrix_rank( D[0]) == order and la.matrix_rank(sumD) == order - 1: alpha = CTMCSolve(sumD) if np.min(np.abs(alpha)) > prec: fullZero = False for Di in D: if np.all(Di == 0.0): fullZero = True break if not fullZero: # scale to the mean value m = MarginalMomentsFromMMAP(D, 1)[0] for i in range(types + 1): D[i] *= m / mean return D trials += 1 raise Exception( "No feasible random MAP/MMAP found with such many zero entries! Try to increase the maxTrials parameter!" )