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

Generated on Mon 24-Jun-2019 15:58:45 by m2html © 2005