Home > matpower5.1 > toggle_dcline.m

toggle_dcline

PURPOSE ^

TOGGLE_DCLINE Enable, disable or check status of DC line modeling.

SYNOPSIS ^

function mpc = toggle_dcline(mpc, on_off)

DESCRIPTION ^

TOGGLE_DCLINE Enable, disable or check status of DC line modeling.
   MPC = TOGGLE_DCLINE(MPC, 'on')
   MPC = TOGGLE_DCLINE(MPC, 'off')
   T_F = TOGGLE_DCLINE(MPC, 'status')

   Enables, disables or checks the status of a set of OPF userfcn
   callbacks to implement DC lines as a pair of linked generators.
   While it uses the OPF extension mechanism, this implementation
   works for simple power flow as well as OPF problems.

   These callbacks expect to find a 'dcline' field in the input MPC,
   where MPC.dcline is an ndc x 17 matrix with columns as defined
   in IDX_DCLINE, where ndc is the number of DC lines.

   The 'int2ext' callback also packages up flow results and stores them
   in appropriate columns of MPC.dcline.

   NOTE: Because of the way this extension modifies the number of
   rows in the gen and gencost matrices, caution must be taken
   when using it with other extensions that deal with generators.

   Examples:
       mpc = loadcase('t_case9_dcline');
       mpc = toggle_dcline(mpc, 'on');
       results1 = runpf(mpc);
       results2 = runopf(mpc);

   See also IDX_DCLINE, ADD_USERFCN, REMOVE_USERFCN, RUN_USERFCN.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SUBFUNCTIONS ^

SOURCE CODE ^

0001 function mpc = toggle_dcline(mpc, on_off)
0002 %TOGGLE_DCLINE Enable, disable or check status of DC line modeling.
0003 %   MPC = TOGGLE_DCLINE(MPC, 'on')
0004 %   MPC = TOGGLE_DCLINE(MPC, 'off')
0005 %   T_F = TOGGLE_DCLINE(MPC, 'status')
0006 %
0007 %   Enables, disables or checks the status of a set of OPF userfcn
0008 %   callbacks to implement DC lines as a pair of linked generators.
0009 %   While it uses the OPF extension mechanism, this implementation
0010 %   works for simple power flow as well as OPF problems.
0011 %
0012 %   These callbacks expect to find a 'dcline' field in the input MPC,
0013 %   where MPC.dcline is an ndc x 17 matrix with columns as defined
0014 %   in IDX_DCLINE, where ndc is the number of DC lines.
0015 %
0016 %   The 'int2ext' callback also packages up flow results and stores them
0017 %   in appropriate columns of MPC.dcline.
0018 %
0019 %   NOTE: Because of the way this extension modifies the number of
0020 %   rows in the gen and gencost matrices, caution must be taken
0021 %   when using it with other extensions that deal with generators.
0022 %
0023 %   Examples:
0024 %       mpc = loadcase('t_case9_dcline');
0025 %       mpc = toggle_dcline(mpc, 'on');
0026 %       results1 = runpf(mpc);
0027 %       results2 = runopf(mpc);
0028 %
0029 %   See also IDX_DCLINE, ADD_USERFCN, REMOVE_USERFCN, RUN_USERFCN.
0030 
0031 %   MATPOWER
0032 %   Copyright (c) 2011-2015 by Power System Engineering Research Center (PSERC)
0033 %   by Ray Zimmerman, PSERC Cornell
0034 %
0035 %   $Id: toggle_dcline.m 2644 2015-03-11 19:34:22Z ray $
0036 %
0037 %   This file is part of MATPOWER.
0038 %   Covered by the 3-clause BSD License (see LICENSE file for details).
0039 %   See http://www.pserc.cornell.edu/matpower/ for more info.
0040 
0041 if strcmp(upper(on_off), 'ON')
0042     %% define named indices into data matrices
0043     c = idx_dcline;
0044 
0045     %% check for proper input data
0046     if ~isfield(mpc, 'dcline') || size(mpc.dcline, 2) < c.LOSS1
0047         error('toggle_dcline: case must contain a ''dcline'' field, an ndc x %d matrix.', c.LOSS1);
0048     end
0049     if isfield(mpc, 'dclinecost') && size(mpc.dcline, 1) ~= size(mpc.dclinecost, 1)
0050         error('toggle_dcline: number of rows in ''dcline'' field (%d) and ''dclinecost'' field (%d) do not match.', ...
0051             size(mpc.dcline, 1), size(mpc.dclinecost, 1));
0052     end
0053 %     k = find(mpc.dcline(:, c.LOSS1) < 0);
0054 %     if ~isempty(k)
0055 %         warning('toggle_dcline: linear loss term is negative for DC line from bus %d to %d\n', ...
0056 %             [mpc.dcline(k, c.F_BUS:c.T_BUS)]');
0057 %     end
0058 
0059     %% add callback functions
0060     %% note: assumes all necessary data included in 1st arg (mpc, om, results)
0061     %%       so, no additional explicit args are needed
0062     mpc = add_userfcn(mpc, 'ext2int', @userfcn_dcline_ext2int);
0063     mpc = add_userfcn(mpc, 'formulation', @userfcn_dcline_formulation);
0064     mpc = add_userfcn(mpc, 'int2ext', @userfcn_dcline_int2ext);
0065     mpc = add_userfcn(mpc, 'printpf', @userfcn_dcline_printpf);
0066     mpc = add_userfcn(mpc, 'savecase', @userfcn_dcline_savecase);
0067     mpc.userfcn.status.dcline = 1;
0068 elseif strcmp(upper(on_off), 'OFF')
0069     mpc = remove_userfcn(mpc, 'savecase', @userfcn_dcline_savecase);
0070     mpc = remove_userfcn(mpc, 'printpf', @userfcn_dcline_printpf);
0071     mpc = remove_userfcn(mpc, 'int2ext', @userfcn_dcline_int2ext);
0072     mpc = remove_userfcn(mpc, 'formulation', @userfcn_dcline_formulation);
0073     mpc = remove_userfcn(mpc, 'ext2int', @userfcn_dcline_ext2int);
0074     mpc.userfcn.status.dcline = 0;
0075 elseif strcmp(upper(on_off), 'STATUS')
0076     if isfield(mpc, 'userfcn') && isfield(mpc.userfcn, 'status') && ...
0077             isfield(mpc.userfcn.status, 'dcline')
0078         mpc = mpc.userfcn.status.dcline;
0079     else
0080         mpc = 0;
0081     end
0082 else
0083     error('toggle_dcline: 2nd argument must be ''on'', ''off'' or ''status''');
0084 end
0085 
0086 
0087 %%-----  ext2int  ------------------------------------------------------
0088 function mpc = userfcn_dcline_ext2int(mpc, args)
0089 %
0090 %   mpc = userfcn_dcline_ext2int(mpc, args)
0091 %
0092 %   This is the 'ext2int' stage userfcn callback that prepares the input
0093 %   data for the formulation stage. It expects to find a 'dcline' field
0094 %   in mpc as described above. The optional args are not currently used.
0095 %   It adds two dummy generators for each in-service DC line, with the
0096 %   appropriate upper and lower generation bounds and corresponding
0097 %   entries in gencost. It also expands columns of any A and N matrices
0098 %   accordingly, if present.
0099 
0100 %% define named indices into data matrices
0101 [PQ, PV, REF, NONE, BUS_I, BUS_TYPE, PD, QD, GS, BS, BUS_AREA, VM, ...
0102     VA, BASE_KV, ZONE, VMAX, VMIN, LAM_P, LAM_Q, MU_VMAX, MU_VMIN] = idx_bus;
0103 [GEN_BUS, PG, QG, QMAX, QMIN, VG, MBASE, GEN_STATUS, PMAX, PMIN, ...
0104     MU_PMAX, MU_PMIN, MU_QMAX, MU_QMIN, PC1, PC2, QC1MIN, QC1MAX, ...
0105     QC2MIN, QC2MAX, RAMP_AGC, RAMP_10, RAMP_30, RAMP_Q, APF] = idx_gen;
0106 [PW_LINEAR, POLYNOMIAL, MODEL, STARTUP, SHUTDOWN, NCOST, COST] = idx_cost;
0107 c = idx_dcline;
0108 
0109 %% initialize some things
0110 if isfield(mpc, 'dclinecost')
0111     havecost = 1;
0112 else
0113     havecost = 0;
0114 end
0115 
0116 %% save version with external indexing
0117 mpc.order.ext.dcline = mpc.dcline;              %% external indexing
0118 if havecost
0119     mpc.order.ext.dclinecost = mpc.dclinecost;  %% external indexing
0120 end
0121 
0122 %% work with only in-service DC lines
0123 mpc.order.dcline.status.on  = find(mpc.dcline(:, c.BR_STATUS) >  0);
0124 mpc.order.dcline.status.off = find(mpc.dcline(:, c.BR_STATUS) <= 0);
0125 
0126 %% remove out-of-service DC lines
0127 dc = mpc.dcline(mpc.order.dcline.status.on, :); %% only in-service DC lines
0128 if havecost
0129     dcc = mpc.dclinecost(mpc.order.dcline.status.on, :);    %% only in-service DC lines
0130     mpc.dclinecost = dcc;
0131 end
0132 ndc = size(dc, 1);          %% number of in-service DC lines
0133 o = mpc.order;
0134 
0135 %%-----  convert stuff to internal indexing  -----
0136 dc(:, c.F_BUS) = o.bus.e2i(dc(:, c.F_BUS));
0137 dc(:, c.T_BUS) = o.bus.e2i(dc(:, c.T_BUS));
0138 mpc.dcline = dc;
0139 
0140 %%-----  create gens to represent DC line terminals  -----
0141 %% ensure consistency of initial values of PF, PT and losses
0142 %% (for simple power flow cases)
0143 dc(:, c.PT) = dc(:, c.PF) - (dc(:, c.LOSS0) + dc(:, c.LOSS1) .* dc(:, c.PF));
0144 
0145 %% create gens
0146 fg = zeros(ndc, size(mpc.gen, 2));
0147 fg(:, MBASE)        = 100;
0148 fg(:, GEN_STATUS)   =  dc(:, c.BR_STATUS);  %% status (should be all 1's)
0149 fg(:, PMIN)         = -Inf;
0150 fg(:, PMAX)         =  Inf;
0151 tg = fg;
0152 fg(:, GEN_BUS)      =  dc(:, c.F_BUS);      %% from bus
0153 tg(:, GEN_BUS)      =  dc(:, c.T_BUS);      %% to bus
0154 fg(:, PG)           = -dc(:, c.PF);         %% flow (extracted at "from")
0155 tg(:, PG)           =  dc(:, c.PT);         %% flow (injected at "to")
0156 fg(:, QG)           =  dc(:, c.QF);         %% VAr injection at "from"
0157 tg(:, QG)           =  dc(:, c.QT);         %% VAr injection at "to"
0158 fg(:, VG)           =  dc(:, c.VF);         %% voltage set-point at "from"
0159 tg(:, VG)           =  dc(:, c.VT);         %% voltage set-point at "to"
0160 k = find(dc(:, c.PMIN) >= 0);           %% min positive direction flow
0161 if ~isempty(k)                              %% contrain at "from" end
0162     fg(k, PMAX)     = -dc(k, c.PMIN);       %% "from" extraction lower lim
0163 end
0164 k = find(dc(:, c.PMAX) >= 0);           %% max positive direction flow
0165 if ~isempty(k)                              %% contrain at "from" end
0166     fg(k, PMIN)     = -dc(k, c.PMAX);       %% "from" extraction upper lim
0167 end
0168 k = find(dc(:, c.PMIN) < 0);            %% max negative direction flow
0169 if ~isempty(k)                              %% contrain at "to" end
0170     tg(k, PMIN)     =  dc(k, c.PMIN);       %% "to" injection lower lim
0171 end
0172 k = find(dc(:, c.PMAX) < 0);            %% min negative direction flow
0173 if ~isempty(k)                              %% contrain at "to" end
0174     tg(k, PMAX)     =  dc(k, c.PMAX);       %% "to" injection upper lim
0175 end
0176 fg(:, QMIN)         =  dc(:, c.QMINF);      %% "from" VAr injection lower lim
0177 fg(:, QMAX)         =  dc(:, c.QMAXF);      %% "from" VAr injection upper lim
0178 tg(:, QMIN)         =  dc(:, c.QMINT);      %%  "to"  VAr injection lower lim
0179 tg(:, QMAX)         =  dc(:, c.QMAXT);      %%  "to"  VAr injection upper lim
0180 
0181 %% fudge PMAX a bit if necessary to avoid triggering
0182 %% dispatchable load constant power factor constraints
0183 fg(isload(fg), PMAX) = -1e-6;
0184 tg(isload(tg), PMAX) = -1e-6;
0185 
0186 %% set all terminal buses to PV (except ref bus)
0187 refbus = find(mpc.bus(:, BUS_TYPE) == REF);
0188 mpc.bus(dc(:, c.F_BUS), BUS_TYPE) = PV;
0189 mpc.bus(dc(:, c.T_BUS), BUS_TYPE) = PV;
0190 mpc.bus(refbus, BUS_TYPE) = REF;
0191 
0192 %% expand A and N, if present
0193 nb = size(mpc.bus, 1);
0194 ng = size(mpc.gen, 1);
0195 if isfield(mpc, 'A') && ~isempty(mpc.A)
0196     [mA, nA] = size(mpc.A);
0197     if nA >= 2*nb + 2*ng    %% assume AC dimensions
0198         mpc.A = [   mpc.A(:, 1:2*nb+ng)      sparse(mA, 2*ndc) ...
0199                     mpc.A(:, 2*nb+ng+(1:ng)) sparse(mA, 2*ndc) ...
0200                     mpc.A(:, (2*nb+2*ng+1:nA)) ];
0201     else                    %% assume DC dimensions
0202         mpc.A = [   mpc.A(:, 1:nb+ng)   sparse(mA, 2*ndc) ...
0203                     mpc.A(:, (nb+ng+1:nA)) ];
0204     end
0205 end
0206 if isfield(mpc, 'N') && ~isempty(mpc.N)
0207     [mN, nN] = size(mpc.N);
0208     if nN >= 2*nb + 2*ng    %% assume AC dimensions
0209         mpc.N = [   mpc.N(:, 1:2*nb+ng)      sparse(mN, 2*ndc) ...
0210                     mpc.N(:, 2*nb+ng+(1:ng)) sparse(mN, 2*ndc) ...
0211                     mpc.N(:, (2*nb+2*ng+1:nN)) ];
0212     else                    %% assume DC dimensions
0213         mpc.N = [   mpc.N(:, 1:nb+ng)   sparse(mN, 2*ndc) ...
0214                     mpc.N(:, (nb+ng+1:nN)) ];
0215     end
0216 end
0217 
0218 %% append dummy gens
0219 mpc.gen = [mpc.gen; fg; tg];
0220 
0221 %% gencost
0222 if isfield(mpc, 'gencost') && ~isempty(mpc.gencost)
0223     [ngcr, ngcc] = size(mpc.gencost);   %% dimensions of gencost
0224     if havecost         %% user has provided costs
0225         ndccc = size(dcc, 2);           %% number of dclinecost columns
0226         ccc = max([ngcc; ndccc]);       %% number of columns in new gencost
0227         if ccc > ngcc                   %% right zero-pad gencost
0228             mpc.gencost = [mpc.gencost zeros(ngcr, ccc-ngcc)];
0229         end
0230 
0231         %% flip function across vertical axis and append to gencost
0232         %% (PF for DC line = -PG for dummy gen at "from" bus)
0233         for k = 1:ndc
0234             if dcc(k, MODEL) == POLYNOMIAL
0235                 nc = dcc(k, NCOST);
0236                 temp = dcc(k, NCOST+(1:nc));
0237                 %% flip sign on coefficients of odd terms
0238                 %% (every other starting with linear term,
0239                 %%  that is, the next to last one)
0240                 temp((nc-1):-2:1) = -temp((nc-1):-2:1);
0241             else  %% dcc(k, MODEL) == PW_LINEAR
0242                 nc = dcc(k, NCOST);
0243                 temp = dcc(k, NCOST+(1:2*nc));
0244                 %% switch sign on horizontal coordinate
0245                 xx = -temp(1:2:2*nc);
0246                 yy =  temp(2:2:2*nc);
0247                 temp(1:2:2*nc) = xx(end:-1:1);
0248                 temp(2:2:2*nc) = yy(end:-1:1);
0249             end
0250             padding = zeros(1, ccc-NCOST-length(temp));
0251             gck = [dcc(k, 1:NCOST) temp padding];
0252             
0253             %% append to gencost
0254             mpc.gencost = [mpc.gencost; gck];
0255         end
0256         %% use zero cost on "to" end gen
0257         tgc = ones(ndc, 1) * [2 0 0 2 zeros(1, ccc-4)];
0258         mpc.gencost = [mpc.gencost; tgc];        
0259     else
0260         %% use zero cost as default
0261         dcgc = ones(2*ndc, 1) * [2 0 0 2 zeros(1, ngcc-4)];
0262         mpc.gencost = [mpc.gencost; dcgc];
0263     end
0264 end
0265 
0266 
0267 %%-----  formulation  --------------------------------------------------
0268 function om = userfcn_dcline_formulation(om, args)
0269 %
0270 %   om = userfcn_dcline_formulation(om, args)
0271 %
0272 %   This is the 'formulation' stage userfcn callback that defines the
0273 %   user constraints for the dummy generators representing DC lines.
0274 %   It expects to find a 'dcline' field in the mpc stored in om, as
0275 %   described above. By the time it is passed to this callback,
0276 %   MPC.dcline should contain only in-service lines and the from and
0277 %   two bus columns should be converted to internal indexing. The
0278 %   optional args are not currently used.
0279 %
0280 %   If Pf, Pt and Ploss are the flow at the "from" end, flow at the
0281 %   "to" end and loss respectively, and L0 and L1 are the linear loss
0282 %   coefficients, the the relationships between them is given by:
0283 %       Pf - Ploss = Pt
0284 %       Ploss = L0 + L1 * Pf
0285 %   If Pgf and Pgt represent the injections of the dummy generators
0286 %   representing the DC line injections into the network, then
0287 %   Pgf = -Pf and Pgt = Pt, and we can combine all of the above to
0288 %   get the following constraint on Pgf ang Pgt:
0289 %       -Pgf - (L0 - L1 * Pgf) = Pgt
0290 %   which can be written:
0291 %       -L0 <= (1 - L1) * Pgf + Pgt <= -L0
0292 
0293 %% define named indices into data matrices
0294 c = idx_dcline;
0295 
0296 %% initialize some things
0297 mpc = get_mpc(om);
0298 dc = mpc.dcline;
0299 ndc = size(dc, 1);              %% number of in-service DC lines
0300 ng  = size(mpc.gen, 1) - 2*ndc; %% number of original gens/disp loads
0301 
0302 %% constraints
0303 nL0 = -dc(:, c.LOSS0) / mpc.baseMVA;
0304 L1  =  dc(:, c.LOSS1);
0305 Adc = [sparse(ndc, ng) spdiags(1-L1, 0, ndc, ndc) speye(ndc, ndc)];
0306 
0307 %% add them to the model
0308 om = add_constraints(om, 'dcline', Adc, nL0, nL0, {'Pg'});
0309 
0310 
0311 %%-----  int2ext  ------------------------------------------------------
0312 function results = userfcn_dcline_int2ext(results, args)
0313 %
0314 %   results = userfcn_dcline_int2ext(results, args)
0315 %
0316 %   This is the 'int2ext' stage userfcn callback that converts everything
0317 %   back to external indexing and packages up the results. It expects to
0318 %   find a 'dcline' field in the results struct as described for mpc
0319 %   above. It also expects that the last 2*ndc entries in the gen and
0320 %   gencost matrices correspond to the in-service DC lines (where ndc is
0321 %   the number of rows in MPC.dcline. These extra rows are removed from
0322 %   gen and gencost and the flow is taken from the PG of these gens and
0323 %   placed in the flow column of the appropiate dcline row. Corresponding
0324 %   columns are also removed from any A and N matrices, if present. The
0325 %   optional args are not currently used.
0326 
0327 %% define named indices into data matrices
0328 [GEN_BUS, PG, QG, QMAX, QMIN, VG, MBASE, GEN_STATUS, PMAX, PMIN, ...
0329     MU_PMAX, MU_PMIN, MU_QMAX, MU_QMIN, PC1, PC2, QC1MIN, QC1MAX, ...
0330     QC2MIN, QC2MAX, RAMP_AGC, RAMP_10, RAMP_30, RAMP_Q, APF] = idx_gen;
0331 c = idx_dcline;
0332 
0333 %% initialize some things
0334 o = results.order;
0335 k = find(o.ext.dcline(:, c.BR_STATUS));
0336 ndc = length(k);                    %% number of in-service DC lines
0337 ng  = size(results.gen, 1) - 2*ndc; %% number of original gens/disp loads
0338 nb  = size(results.bus, 1);
0339 
0340 %% extract dummy gens
0341 fg = results.gen(ng    +(1:ndc), :);
0342 tg = results.gen(ng+ndc+(1:ndc), :);
0343 
0344 %% remove dummy gens
0345 results.gen     = results.gen(1:ng, :);
0346 if isfield(results, 'gencost') && ~isempty(results.gencost)
0347     results.gencost = results.gencost(1:ng, :);
0348 end
0349 
0350 %% delete corresponding columns from A and N, if present
0351 if isfield(results, 'A') && ~isempty(results.A)
0352     [mA, nA] = size(results.A);
0353     if nA >= 2*nb + 2*ng + 4*ndc    %% assume AC dimensions
0354         results.A = results.A(:, [1:2*nb+ng 2*nb+ng+2*ndc+(1:ng) 2*nb+2*ng+4*ndc+1:nA]);
0355     else                            %% assume DC dimensions
0356         results.A = results.A(:, [1:nb+ng nb+ng+2*ndc+1:nA]);
0357     end
0358 end
0359 if isfield(results, 'N') && ~isempty(results.N)
0360     [mN, nN] = size(results.N);
0361     if nN >= 2*nb + 2*ng + 4*ndc    %% assume AC dimensions
0362         results.N = results.N(:, [1:2*nb+ng 2*nb+ng+2*ndc+(1:ng) 2*nb+2*ng+4*ndc+1:nN]);
0363     else                            %% assume DC dimensions
0364         results.N = results.N(:, [1:nb+ng nb+ng+2*ndc+1:nN]);
0365     end
0366 end
0367 
0368 %% get the solved flows
0369 results.dcline(:, c.PF) = -fg(:, PG);
0370 results.dcline(:, c.PT) =  tg(:, PG);
0371 results.dcline(:, c.QF) =  fg(:, QG);
0372 results.dcline(:, c.QT) =  tg(:, QG);
0373 results.dcline(:, c.VF) =  fg(:, VG);
0374 results.dcline(:, c.VT) =  tg(:, VG);
0375 if size(fg, 2) >= MU_QMIN
0376     results.dcline(:, c.MU_PMIN ) = fg(:, MU_PMAX) + tg(:, MU_PMIN);
0377     results.dcline(:, c.MU_PMAX ) = fg(:, MU_PMIN) + tg(:, MU_PMAX);
0378     results.dcline(:, c.MU_QMINF) = fg(:, MU_QMIN);
0379     results.dcline(:, c.MU_QMAXF) = fg(:, MU_QMAX);
0380     results.dcline(:, c.MU_QMINT) = tg(:, MU_QMIN);
0381     results.dcline(:, c.MU_QMAXT) = tg(:, MU_QMAX);
0382 end
0383 
0384 %%-----  convert stuff back to external indexing  -----
0385 results.order.int.dcline = results.dcline;  %% save internal version
0386 %% copy results to external version
0387 o.ext.dcline(k, c.PF:c.VT) = results.dcline(:, c.PF:c.VT);
0388 if size(results.dcline, 2) == c.MU_QMAXT
0389     o.ext.dcline(k, c.MU_PMIN:c.MU_QMAXT) = results.dcline(:, c.MU_PMIN:c.MU_QMAXT);
0390 end
0391 results.dcline = o.ext.dcline;              %% use external version
0392 
0393 
0394 %%-----  printpf  ------------------------------------------------------
0395 function results = userfcn_dcline_printpf(results, fd, mpopt, args)
0396 %
0397 %   results = userfcn_dcline_printpf(results, fd, mpopt, args)
0398 %
0399 %   This is the 'printpf' stage userfcn callback that pretty-prints the
0400 %   results. It expects a results struct, a file descriptor and a MATPOWER
0401 %   options struct. The optional args are not currently used.
0402 
0403 %% define named indices into data matrices
0404 c = idx_dcline;
0405 
0406 %% options
0407 SUPPRESS        = mpopt.out.suppress_detail;
0408 if SUPPRESS == -1
0409     if size(results.bus, 1) > 500
0410         SUPPRESS = 1;
0411     else
0412         SUPPRESS = 0;
0413     end
0414 end
0415 OUT_ALL = mpopt.out.all;
0416 OUT_BRANCH      = OUT_ALL == 1 || (OUT_ALL == -1 && ~SUPPRESS && mpopt.out.branch);
0417 if OUT_ALL == -1
0418     OUT_ALL_LIM = ~SUPPRESS * mpopt.out.lim.all;
0419 elseif OUT_ALL == 1
0420     OUT_ALL_LIM = 2;
0421 else
0422     OUT_ALL_LIM = 0;
0423 end
0424 if OUT_ALL_LIM == -1
0425     OUT_LINE_LIM    = ~SUPPRESS * mpopt.out.lim.line;
0426 else
0427     OUT_LINE_LIM    = OUT_ALL_LIM;
0428 end
0429 ctol = mpopt.opf.violation; %% constraint violation tolerance
0430 ptol = 1e-4;                %% tolerance for displaying shadow prices
0431 
0432 %%-----  print results  -----
0433 dc = results.dcline;
0434 ndc = size(dc, 1);
0435 kk = find(dc(:, c.BR_STATUS) ~= 0);
0436 if OUT_BRANCH
0437     fprintf(fd, '\n================================================================================');
0438     fprintf(fd, '\n|     DC Line Data                                                             |');
0439     fprintf(fd, '\n================================================================================');
0440     fprintf(fd, '\n Line    From     To        Power Flow           Loss     Reactive Inj (MVAr)');
0441     fprintf(fd, '\n   #      Bus     Bus   From (MW)   To (MW)      (MW)       From        To   ');
0442     fprintf(fd, '\n------  ------  ------  ---------  ---------  ---------  ---------  ---------');
0443     loss = 0;
0444     for k = 1:ndc
0445         if dc(k, c.BR_STATUS)   %% status on
0446             fprintf(fd, '\n%5d%8d%8d%11.2f%11.2f%11.2f%11.2f%11.2f', ...
0447                         k, dc(k, c.F_BUS:c.T_BUS), dc(k, c.PF:c.PT), ...
0448                         dc(k, c.PF) - dc(k, c.PT), dc(k, c.QF:c.QT) );
0449             loss = loss + dc(k, c.PF) - dc(k, c.PT);
0450         else
0451             fprintf(fd, '\n%5d%8d%8d%11s%11s%11s%11s%11s', ...
0452                         k, dc(k, c.F_BUS:c.T_BUS), '-  ', '-  ', '-  ', '-  ', '-  ');
0453         end
0454     end
0455     fprintf(fd, '\n                                              ---------');
0456     fprintf(fd, '\n                                     Total:%11.2f\n', loss);
0457 end
0458 
0459 if OUT_LINE_LIM == 2 || (OUT_LINE_LIM == 1 && ...
0460         (any(dc(kk, c.PF) > dc(kk, c.PMAX) - ctol) || ...
0461          any(dc(kk, c.MU_PMIN) > ptol) || ...
0462          any(dc(kk, c.MU_PMAX) > ptol)))
0463     fprintf(fd, '\n================================================================================');
0464     fprintf(fd, '\n|     DC Line Constraints                                                      |');
0465     fprintf(fd, '\n================================================================================');
0466     fprintf(fd, '\n Line    From     To          Minimum        Actual Flow       Maximum');
0467     fprintf(fd, '\n   #      Bus     Bus    Pmin mu     Pmin       (MW)       Pmax      Pmax mu ');
0468     fprintf(fd, '\n------  ------  ------  ---------  ---------  ---------  ---------  ---------');
0469     for k = 1:ndc
0470         if OUT_LINE_LIM == 2 || (OUT_LINE_LIM == 1 && ...
0471                 (dc(k, c.PF) > dc(k, c.PMAX) - ctol || ...
0472                  dc(k, c.MU_PMIN) > ptol || ...
0473                  dc(k, c.MU_PMAX) > ptol))
0474             if dc(k, c.BR_STATUS)   %% status on
0475                 fprintf(fd, '\n%5d%8d%8d', k, dc(k, c.F_BUS:c.T_BUS) );
0476                 if dc(k, c.MU_PMIN) > ptol
0477                     fprintf(fd, '%11.3f', dc(k, c.MU_PMIN) );
0478                 else
0479                     fprintf(fd, '%11s', '-  ' );
0480                 end
0481                 fprintf(fd, '%11.2f%11.2f%11.2f', ...
0482                             dc(k, c.PMIN), dc(k, c.PF), dc(k, c.PMAX) );
0483                 if dc(k, c.MU_PMAX) > ptol
0484                     fprintf(fd, '%11.3f', dc(k, c.MU_PMAX) );
0485                 else
0486                     fprintf(fd, '%11s', '-  ' );
0487                 end
0488             else
0489                 fprintf(fd, '\n%5d%8d%8d%11s%11s%11s%11s%11s', ...
0490                             k, dc(k, c.F_BUS:c.T_BUS), '-  ', '-  ', '-  ', '-  ', '-  ');
0491             end
0492         end
0493     end
0494     fprintf(fd, '\n');
0495 end
0496 
0497 
0498 %%-----  savecase  -----------------------------------------------------
0499 function mpc = userfcn_dcline_savecase(mpc, fd, prefix, args)
0500 %
0501 %   mpc = userfcn_dcline_savecase(mpc, fd, mpopt, args)
0502 %
0503 %   This is the 'savecase' stage userfcn callback that prints the M-file
0504 %   code to save the 'dcline' field in the case file. It expects a
0505 %   MATPOWER case struct (mpc), a file descriptor and variable prefix
0506 %   (usually 'mpc.'). The optional args are not currently used.
0507 
0508 %% define named indices into data matrices
0509 c = idx_dcline;
0510 
0511 %% save it
0512 ncols = size(mpc.dcline, 2);
0513 fprintf(fd, '\n%%%%-----  DC Line Data  -----%%%%\n');
0514 if ncols < c.MU_QMAXT
0515     fprintf(fd, '%%\tfbus\ttbus\tstatus\tPf\tPt\tQf\tQt\tVf\tVt\tPmin\tPmax\tQminF\tQmaxF\tQminT\tQmaxT\tloss0\tloss1\n');
0516 else
0517     fprintf(fd, '%%\tfbus\ttbus\tstatus\tPf\tPt\tQf\tQt\tVf\tVt\tPmin\tPmax\tQminF\tQmaxF\tQminT\tQmaxT\tloss0\tloss1\tmuPmin\tmuPmax\tmuQminF\tmuQmaxF\tmuQminT\tmuQmaxT\n');
0518 end
0519 template = '\t%d\t%d\t%d\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g\t%.9g';
0520 if ncols == c.MU_QMAXT
0521     template = [template, '\t%.4f\t%.4f\t%.4f\t%.4f\t%.4f\t%.4f'];
0522 end
0523 template = [template, ';\n'];
0524 fprintf(fd, '%sdcline = [\n', prefix);
0525 fprintf(fd, template, mpc.dcline.');
0526 fprintf(fd, '];\n');
0527 
0528 %% to do, make it save mpc.dclinecost too (somehow I forgot this)
0529 %% have a look at the saving of gencost in savecase for template

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