1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17 """Solves the power flow using a fast decoupled method.
18 """
19
20 import sys
21
22 from numpy import array, angle, exp, linalg, conj, r_, Inf
23 from scipy.sparse.linalg import splu
24
25 from ppoption import ppoption
26
27
28 -def fdpf(Ybus, Sbus, V0, Bp, Bpp, ref, pv, pq, ppopt=None):
29 """Solves the power flow using a fast decoupled method.
30
31 Solves for bus voltages given the full system admittance matrix (for
32 all buses), the complex bus power injection vector (for all buses),
33 the initial vector of complex bus voltages, the FDPF matrices B prime
34 and B double prime, and column vectors with the lists of bus indices
35 for the swing bus, PV buses, and PQ buses, respectively. The bus voltage
36 vector contains the set point for generator (including ref bus)
37 buses, and the reference angle of the swing bus, as well as an initial
38 guess for remaining magnitudes and angles. C{ppopt} is a PYPOWER options
39 vector which can be used to set the termination tolerance, maximum
40 number of iterations, and output options (see L{ppoption} for details).
41 Uses default options if this parameter is not given. Returns the
42 final complex voltages, a flag which indicates whether it converged
43 or not, and the number of iterations performed.
44
45 @see: L{runpf}
46
47 @author: Ray Zimmerman (PSERC Cornell)
48 @author: Richard Lincoln
49 """
50 if ppopt is None:
51 ppopt = ppoption()
52
53
54 tol = ppopt['PF_TOL']
55 max_it = ppopt['PF_MAX_IT_FD']
56 verbose = ppopt['VERBOSE']
57
58
59 converged = 0
60 i = 0
61 V = V0
62 Va = angle(V)
63 Vm = abs(V)
64
65
66
67
68 pvpq = r_[pv, pq]
69
70
71 mis = (V * conj(Ybus * V) - Sbus) / Vm
72 P = mis[pvpq].real
73 Q = mis[pq].imag
74
75
76 normP = linalg.norm(P, Inf)
77 normQ = linalg.norm(Q, Inf)
78 if verbose > 0:
79 alg = ppopt['PF_ALG']
80 s = 'XB' if alg == 2 else 'BX'
81 sys.stdout.write('(fast-decoupled, %s)\n' % s)
82 if verbose > 1:
83 sys.stdout.write('\niteration max mismatch (p.u.) ')
84 sys.stdout.write('\ntype # P Q ')
85 sys.stdout.write('\n---- ---- ----------- -----------')
86 sys.stdout.write('\n - %3d %10.3e %10.3e' % (i, normP, normQ))
87 if normP < tol and normQ < tol:
88 converged = 1
89 if verbose > 1:
90 sys.stdout.write('\nConverged!\n')
91
92
93 Bp = Bp[array([pvpq]).T, pvpq].tocsc()
94 Bpp = Bpp[array([pq]).T, pq].tocsc()
95
96
97 Bp_solver = splu(Bp)
98 Bpp_solver = splu(Bpp)
99
100
101 while (not converged and i < max_it):
102
103 i = i + 1
104
105
106 dVa = -Bp_solver.solve(P)
107
108
109 Va[pvpq] = Va[pvpq] + dVa
110 V = Vm * exp(1j * Va)
111
112
113 mis = (V * conj(Ybus * V) - Sbus) / Vm
114 P = mis[pvpq].real
115 Q = mis[pq].imag
116
117
118 normP = linalg.norm(P, Inf)
119 normQ = linalg.norm(Q, Inf)
120 if verbose > 1:
121 sys.stdout.write("\n %s %3d %10.3e %10.3e" %
122 (type,i, normP, normQ))
123 if normP < tol and normQ < tol:
124 converged = 1
125 if verbose:
126 sys.stdout.write('\nFast-decoupled power flow converged in %d '
127 'P-iterations and %d Q-iterations.\n' % (i, i - 1))
128 break
129
130
131 dVm = -Bpp_solver.solve(Q)
132
133
134 Vm[pq] = Vm[pq] + dVm
135 V = Vm * exp(1j * Va)
136
137
138 mis = (V * conj(Ybus * V) - Sbus) / Vm
139 P = mis[pvpq].real
140 Q = mis[pq].imag
141
142
143 normP = linalg.norm(P, Inf)
144 normQ = linalg.norm(Q, Inf)
145 if verbose > 1:
146 sys.stdout.write('\n Q %3d %10.3e %10.3e' % (i, normP, normQ))
147 if normP < tol and normQ < tol:
148 converged = 1
149 if verbose:
150 sys.stdout.write('\nFast-decoupled power flow converged in %d '
151 'P-iterations and %d Q-iterations.\n' % (i, i))
152 break
153
154 if verbose:
155 if not converged:
156 sys.stdout.write('\nFast-decoupled power flow did not converge in '
157 '%d iterations.' % i)
158
159 return V, converged, i
160