Package pypower :: Module makeYbus
[hide private]
[frames] | no frames]

Source Code for Module pypower.makeYbus

  1  # Copyright (C) 1996-2011 Power System Engineering Research Center (PSERC) 
  2  # Copyright (C) 2011 Richard Lincoln 
  3  # 
  4  # PYPOWER is free software: you can redistribute it and/or modify 
  5  # it under the terms of the GNU General Public License as published 
  6  # by the Free Software Foundation, either version 3 of the License, 
  7  # or (at your option) any later version. 
  8  # 
  9  # PYPOWER is distributed in the hope that it will be useful, 
 10  # but WITHOUT ANY WARRANTY; without even the implied warranty of 
 11  # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 
 12  # GNU General Public License for more details. 
 13  # 
 14  # You should have received a copy of the GNU General Public License 
 15  # along with PYPOWER. If not, see <http://www.gnu.org/licenses/>. 
 16   
 17  """Builds the bus admittance matrix and branch admittance matrices. 
 18  """ 
 19   
 20  from sys import stderr 
 21   
 22  from numpy import ones, conj, nonzero, any, exp, pi, r_ 
 23  from scipy.sparse import csr_matrix 
 24   
 25  from idx_bus import BUS_I, GS, BS 
 26  from idx_brch import F_BUS, T_BUS, BR_R, BR_X, BR_B, BR_STATUS, SHIFT, TAP 
 27   
 28   
29 -def makeYbus(baseMVA, bus, branch):
30 """Builds the bus admittance matrix and branch admittance matrices. 31 32 Returns the full bus admittance matrix (i.e. for all buses) and the 33 matrices C{Yf} and C{Yt} which, when multiplied by a complex voltage 34 vector, yield the vector currents injected into each line from the 35 "from" and "to" buses respectively of each line. Does appropriate 36 conversions to p.u. 37 38 @see: L{makeSbus} 39 40 @author: Ray Zimmerman (PSERC Cornell) 41 @author: Richard Lincoln 42 """ 43 ## constants 44 nb = bus.shape[0] ## number of buses 45 nl = branch.shape[0] ## number of lines 46 47 ## check that bus numbers are equal to indices to bus (one set of bus nums) 48 if any(bus[:, BUS_I] != range(nb)): 49 stderr.write('buses must appear in order by bus number\n') 50 51 ## for each branch, compute the elements of the branch admittance matrix where 52 ## 53 ## | If | | Yff Yft | | Vf | 54 ## | | = | | * | | 55 ## | It | | Ytf Ytt | | Vt | 56 ## 57 stat = branch[:, BR_STATUS] ## ones at in-service branches 58 Ys = stat / (branch[:, BR_R] + 1j * branch[:, BR_X]) ## series admittance 59 Bc = stat * branch[:, BR_B] ## line charging susceptance 60 tap = ones(nl) ## default tap ratio = 1 61 i = nonzero(branch[:, TAP]) ## indices of non-zero tap ratios 62 tap[i] = branch[i, TAP] ## assign non-zero tap ratios 63 tap = tap * exp(1j * pi / 180 * branch[:, SHIFT]) ## add phase shifters 64 65 Ytt = Ys + 1j * Bc / 2 66 Yff = Ytt / (tap * conj(tap)) 67 Yft = - Ys / conj(tap) 68 Ytf = - Ys / tap 69 70 ## compute shunt admittance 71 ## if Psh is the real power consumed by the shunt at V = 1.0 p.u. 72 ## and Qsh is the reactive power injected by the shunt at V = 1.0 p.u. 73 ## then Psh - j Qsh = V * conj(Ysh * V) = conj(Ysh) = Gs - j Bs, 74 ## i.e. Ysh = Psh + j Qsh, so ... 75 ## vector of shunt admittances 76 Ysh = (bus[:, GS] + 1j * bus[:, BS]) / baseMVA 77 78 ## build connection matrices 79 f = branch[:, F_BUS] ## list of "from" buses 80 t = branch[:, T_BUS] ## list of "to" buses 81 ## connection matrix for line & from buses 82 Cf = csr_matrix((ones(nl), (range(nl), f)), (nl, nb)) 83 ## connection matrix for line & to buses 84 Ct = csr_matrix((ones(nl), (range(nl), t)), (nl, nb)) 85 86 ## build Yf and Yt such that Yf * V is the vector of complex branch currents injected 87 ## at each branch's "from" bus, and Yt is the same for the "to" bus end 88 i = r_[range(nl), range(nl)] ## double set of row indices 89 90 Yf = csr_matrix((r_[Yff, Yft], (i, r_[f, t])), (nl, nb)) 91 Yt = csr_matrix((r_[Ytf, Ytt], (i, r_[f, t])), (nl, nb)) 92 # Yf = spdiags(Yff, 0, nl, nl) * Cf + spdiags(Yft, 0, nl, nl) * Ct 93 # Yt = spdiags(Ytf, 0, nl, nl) * Cf + spdiags(Ytt, 0, nl, nl) * Ct 94 95 ## build Ybus 96 Ybus = Cf.T * Yf + Ct.T * Yt + \ 97 csr_matrix((Ysh, (range(nb), range(nb))), (nb, nb)) 98 99 return Ybus, Yf, Yt
100