Home > matpower5.1 > extras > state_estimator > runse.m

runse

PURPOSE ^

RUNSE Runs a state estimator.

SYNOPSIS ^

function [MVAbase, bus, gen, branch, success, et] = runse(casedata, mpopt, fname, solvedcase)

DESCRIPTION ^

RUNSE  Runs a state estimator.
   [BASEMVA, BUS, GEN, BRANCH, SUCCESS, ET] = ...
           RUNSE(CASEDATA, MPOPT, FNAME, SOLVEDCASE)

   Runs a state estimator (after a Newton power flow). Under construction with
   parts based on code from James S. Thorp.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [MVAbase, bus, gen, branch, success, et] = runse(casedata, mpopt, fname, solvedcase)
0002 %RUNSE  Runs a state estimator.
0003 %   [BASEMVA, BUS, GEN, BRANCH, SUCCESS, ET] = ...
0004 %           RUNSE(CASEDATA, MPOPT, FNAME, SOLVEDCASE)
0005 %
0006 %   Runs a state estimator (after a Newton power flow). Under construction with
0007 %   parts based on code from James S. Thorp.
0008 
0009 %   MATPOWER
0010 %   Copyright (c) 1996-2015 by Power System Engineering Research Center (PSERC)
0011 %   by Ray Zimmerman, PSERC Cornell
0012 %   parts based on code by James S. Thorp, June 2004
0013 %
0014 %   $Id: runse.m 2644 2015-03-11 19:34:22Z ray $
0015 %
0016 %   This file is part of MATPOWER.
0017 %   Covered by the 3-clause BSD License (see LICENSE file for details).
0018 %   See http://www.pserc.cornell.edu/matpower/ for more info.
0019 
0020 %%-----  initialize  -----
0021 %% define named indices into bus, gen, branch matrices
0022 [PQ, PV, REF, NONE, BUS_I, BUS_TYPE, PD, QD, GS, BS, BUS_AREA, VM, ...
0023     VA, BASE_KV, ZONE, VMAX, VMIN, LAM_P, LAM_Q, MU_VMAX, MU_VMIN] = idx_bus;
0024 [F_BUS, T_BUS, BR_R, BR_X, BR_B, RATE_A, RATE_B, RATE_C, ...
0025     TAP, SHIFT, BR_STATUS, PF, QF, PT, QT, MU_SF, MU_ST, ...
0026     ANGMIN, ANGMAX, MU_ANGMIN, MU_ANGMAX] = idx_brch;
0027 [GEN_BUS, PG, QG, QMAX, QMIN, VG, MBASE, GEN_STATUS, PMAX, PMIN, ...
0028     MU_PMAX, MU_PMIN, MU_QMAX, MU_QMIN, PC1, PC2, QC1MIN, QC1MAX, ...
0029     QC2MIN, QC2MAX, RAMP_AGC, RAMP_10, RAMP_30, RAMP_Q, APF] = idx_gen;
0030 
0031 %% default arguments
0032 if nargin < 4
0033     solvedcase = '';                %% don't save solved case
0034     if nargin < 3
0035         fname = '';                 %% don't print results to a file
0036         if nargin < 2
0037             mpopt = mpoption;       %% use default options
0038             if nargin < 1
0039                 casedata = 'case9'; %% default data file is 'case9.m'
0040             end
0041         end
0042     end
0043 end
0044 
0045 %% options
0046 dc = strcmp(upper(mpopt.model), 'DC');  %% use DC formulation?
0047 
0048 %% read data & convert to internal bus numbering
0049 [baseMVA, bus, gen, branch] = loadcase(casedata);
0050 [i2e, bus, gen, branch] = ext2int(bus, gen, branch);
0051 
0052 %% get bus index lists of each type of bus
0053 [ref, pv, pq] = bustypes(bus, gen);
0054 
0055 %% generator info
0056 on = find(gen(:, GEN_STATUS) > 0);      %% which generators are on?
0057 gbus = gen(on, GEN_BUS);                %% what buses are they at?
0058 
0059 %%-----  run the power flow  -----
0060 t0 = clock;
0061 if dc                               %% DC formulation
0062     %% initial state
0063     Va0 = bus(:, VA) * (pi/180);
0064     
0065     %% build B matrices and phase shift injections
0066     [B, Bf, Pbusinj, Pfinj] = makeBdc(baseMVA, bus, branch);
0067     
0068     %% compute complex bus power injections (generation - load)
0069     %% adjusted for phase shifters and real shunts
0070     Pbus = real(makeSbus(baseMVA, bus, gen)) - Pbusinj - bus(:, GS) / baseMVA;
0071     
0072     %% "run" the power flow
0073     Va = dcpf(B, Pbus, Va0, ref, pv, pq);
0074     
0075     %% update data matrices with solution
0076     branch(:, [QF, QT]) = zeros(size(branch, 1), 2);
0077     branch(:, PF) = (Bf * Va + Pfinj) * baseMVA;
0078     branch(:, PT) = -branch(:, PF);
0079     bus(:, VM) = ones(size(bus, 1), 1);
0080     bus(:, VA) = Va * (180/pi);
0081     %% update Pg for swing generator (note: other gens at ref bus are accounted for in Pbus)
0082     %%      Pg = Pinj + Pload + Gs
0083     %%      newPg = oldPg + newPinj - oldPinj
0084     refgen = find(gbus == ref);             %% which is(are) the reference gen(s)?
0085     gen(on(refgen(1)), PG) = gen(on(refgen(1)), PG) + (B(ref, :) * Va - Pbus(ref)) * baseMVA;
0086     
0087     success = 1;
0088 else                                %% AC formulation
0089     %% initial state
0090     % V0    = ones(size(bus, 1), 1);            %% flat start
0091     V0  = bus(:, VM) .* exp(sqrt(-1) * pi/180 * bus(:, VA));
0092     V0(gbus) = gen(on, VG) ./ abs(V0(gbus)).* V0(gbus);
0093     
0094     %% build admittance matrices
0095     [Ybus, Yf, Yt] = makeYbus(baseMVA, bus, branch);
0096     
0097     %% compute complex bus power injections (generation - load)
0098     Sbus = makeSbus(baseMVA, bus, gen);
0099     
0100     %% run the power flow
0101     alg = upper(mpopt.pf.alg);
0102     switch alg
0103         case 'NR'
0104             [V, success, iterations] = newtonpf(Ybus, Sbus, V0, ref, pv, pq, mpopt);
0105         case {'FDXB', 'FDBX'}
0106             [Bp, Bpp] = makeB(baseMVA, bus, branch, alg);
0107             [V, success, iterations] = fdpf(Ybus, Sbus, V0, Bp, Bpp, ref, pv, pq, mpopt);
0108         case 'GS'
0109             [V, success, iterations] = gausspf(Ybus, Sbus, V0, ref, pv, pq, mpopt);
0110         otherwise
0111             error('Only Newton''s method, fast-decoupled, and Gauss-Seidel power flow algorithms currently implemented.');
0112     end
0113     
0114     %% update data matrices with solution
0115     [bus, gen, branch] = pfsoln(baseMVA, bus, gen, branch, Ybus, Yf, Yt, V, ref, pv, pq);
0116 end
0117 et = etime(clock, t0);
0118 
0119 %%--------------------  begin state estimator code  --------------------
0120 %% save some values from load flow solution
0121 Pflf=branch(:,PF);
0122 Qflf=branch(:,QF);
0123 Ptlf=branch(:,PT);
0124 Qtlf=branch(:,QT);
0125 Sbuslf = V .* conj(Ybus * V);
0126 Vlf=V;
0127 
0128 %% run state estimator
0129 [V, converged, i] = state_est(branch, Ybus, Yf, Yt, Sbuslf, Vlf, ref, pv, pq, mpopt);
0130 
0131 %% update data matrices to match estimator solution ...
0132 %% ... bus injections at PQ buses
0133 Sbus = V .* conj(Ybus * V);
0134 bus(pq, PD) = -real(Sbus(pq)) * baseMVA;
0135 bus(pq, QD) = -imag(Sbus(pq)) * baseMVA;
0136 %% ... gen outputs at PV buses
0137 on = find(gen(:, GEN_STATUS) > 0);      %% which generators are on?
0138 gbus = gen(on, GEN_BUS);                %% what buses are they at?
0139 gen(on, PG) = real(Sbus(gbus)) * baseMVA + bus(gbus, PD);   %% inj P + local Pd
0140 %% ... line flows, reference bus injections, etc.
0141 [bus, gen, branch] = pfsoln(baseMVA, bus, gen, branch, Ybus, Yf, Yt, V, ref, pv, pq);
0142 
0143 %% plot differences from load flow solution
0144 Pfe=branch(:,PF);
0145 Qfe=branch(:,QF);
0146 Pte=branch(:,PT);
0147 Qte=branch(:,QT);
0148 nbr = length(Pfe);
0149 subplot(3,2,1), plot(180/pi*(angle(Vlf)-angle(V)),'.'), title('Voltage Angle (deg)');
0150 subplot(3,2,2), plot(abs(Vlf)-abs(V),'.'), title('Voltage Magnitude (p.u.)');
0151 subplot(3,2,3), plot((1:nbr),(Pfe-Pflf),'r.',(1:nbr),(Pte-Ptlf),'b.'), title('Real Flow (MW)');
0152 subplot(3,2,4), plot((1:nbr),(Qfe-Qflf),'r.',(1:nbr),(Qte-Qtlf),'b.'), title('Reactive Flow (MVAr)');
0153 subplot(3,2,5), plot(baseMVA*real(Sbuslf-Sbus), '.'), title('Real Injection (MW)');
0154 subplot(3,2,6), plot(baseMVA*imag(Sbuslf-Sbus), '.'), title('Reactive Injection (MVAr)');
0155 %%--------------------  end state estimator code  --------------------
0156 
0157 %%-----  output results  -----
0158 %% convert back to original bus numbering & print results
0159 [bus, gen, branch] = int2ext(i2e, bus, gen, branch);
0160 if fname
0161     [fd, msg] = fopen(fname, 'at');
0162     if fd == -1
0163         error(msg);
0164     else
0165         if mpopt.out.all == 0
0166             printpf(baseMVA, bus, gen, branch, [], success, et, fd, ...
0167                 mpoption(mpopt, 'out.all', -1));
0168         else
0169             printpf(baseMVA, bus, gen, branch, [], success, et, fd, mpopt);
0170         end
0171         fclose(fd);
0172     end
0173 end
0174 printpf(baseMVA, bus, gen, branch, [], success, et, 1, mpopt);
0175 
0176 %% save solved case
0177 if solvedcase
0178     savecase(solvedcase, baseMVA, bus, gen, branch);
0179 end
0180 
0181 %% this is just to prevent it from printing baseMVA
0182 %% when called with no output arguments
0183 if nargout, MVAbase = baseMVA; end

Generated on Fri 20-Mar-2015 18:23:34 by m2html © 2005