Home > matpower5.1 > opf_costfcn.m

opf_costfcn

PURPOSE ^

OPF_COSTFCN Evaluates objective function, gradient and Hessian for OPF.

SYNOPSIS ^

function [f, df, d2f] = opf_costfcn(x, om, varargin)

DESCRIPTION ^

OPF_COSTFCN  Evaluates objective function, gradient and Hessian for OPF.
   [F, DF, D2F] = OPF_COSTFCN(X, OM)

   Objective function evaluation routine for AC optimal power flow,
   suitable for use with MIPS or FMINCON. Computes objective function value,
   gradient and Hessian.

   Inputs:
     X : optimization vector
     OM : OPF model object

   Outputs:
     F   : value of objective function
     DF  : (optional) gradient of objective function (column vector)
     D2F : (optional) Hessian of objective function (sparse matrix)

   Examples:
       f = opf_costfcn(x, om);
       [f, df] = opf_costfcn(x, om);
       [f, df, d2f] = opf_costfcn(x, om);

   See also OPF_CONSFCN, OPF_HESSFCN.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [f, df, d2f] = opf_costfcn(x, om, varargin)
0002 %OPF_COSTFCN  Evaluates objective function, gradient and Hessian for OPF.
0003 %   [F, DF, D2F] = OPF_COSTFCN(X, OM)
0004 %
0005 %   Objective function evaluation routine for AC optimal power flow,
0006 %   suitable for use with MIPS or FMINCON. Computes objective function value,
0007 %   gradient and Hessian.
0008 %
0009 %   Inputs:
0010 %     X : optimization vector
0011 %     OM : OPF model object
0012 %
0013 %   Outputs:
0014 %     F   : value of objective function
0015 %     DF  : (optional) gradient of objective function (column vector)
0016 %     D2F : (optional) Hessian of objective function (sparse matrix)
0017 %
0018 %   Examples:
0019 %       f = opf_costfcn(x, om);
0020 %       [f, df] = opf_costfcn(x, om);
0021 %       [f, df, d2f] = opf_costfcn(x, om);
0022 %
0023 %   See also OPF_CONSFCN, OPF_HESSFCN.
0024 
0025 %   MATPOWER
0026 %   Copyright (c) 1996-2015 by Power System Engineering Research Center (PSERC)
0027 %   by Carlos E. Murillo-Sanchez, PSERC Cornell & Universidad Autonoma de Manizales
0028 %   and Ray Zimmerman, PSERC Cornell
0029 %
0030 %   $Id: opf_costfcn.m 2644 2015-03-11 19:34:22Z ray $
0031 %
0032 %   This file is part of MATPOWER.
0033 %   Covered by the 3-clause BSD License (see LICENSE file for details).
0034 %   See http://www.pserc.cornell.edu/matpower/ for more info.
0035 
0036 %%----- initialize -----
0037 %% define named indices into data matrices
0038 [GEN_BUS, PG, QG, QMAX, QMIN, VG, MBASE, GEN_STATUS, PMAX, PMIN, ...
0039     MU_PMAX, MU_PMIN, MU_QMAX, MU_QMIN, PC1, PC2, QC1MIN, QC1MAX, ...
0040     QC2MIN, QC2MAX, RAMP_AGC, RAMP_10, RAMP_30, RAMP_Q, APF] = idx_gen;
0041 [F_BUS, T_BUS, BR_R, BR_X, BR_B, RATE_A, RATE_B, RATE_C, ...
0042     TAP, SHIFT, BR_STATUS, PF, QF, PT, QT, MU_SF, MU_ST, ...
0043     ANGMIN, ANGMAX, MU_ANGMIN, MU_ANGMAX] = idx_brch;
0044 [PW_LINEAR, POLYNOMIAL, MODEL, STARTUP, SHUTDOWN, NCOST, COST] = idx_cost;
0045 
0046 %% unpack data
0047 mpc = get_mpc(om);
0048 [baseMVA, gen, gencost] = deal(mpc.baseMVA, mpc.gen, mpc.gencost);
0049 cp = get_cost_params(om);
0050 [N, Cw, H, dd, rh, kk, mm] = deal(cp.N, cp.Cw, cp.H, cp.dd, ...
0051                                     cp.rh, cp.kk, cp.mm);
0052 vv = get_idx(om);
0053 
0054 %% problem dimensions
0055 ng = size(gen, 1);          %% number of dispatchable injections
0056 ny = getN(om, 'var', 'y');  %% number of piece-wise linear costs
0057 nxyz = length(x);           %% total number of control vars of all types
0058 
0059 %% grab Pg & Qg
0060 Pg = x(vv.i1.Pg:vv.iN.Pg);  %% active generation in p.u.
0061 Qg = x(vv.i1.Qg:vv.iN.Qg);  %% reactive generation in p.u.
0062 
0063 %%----- evaluate objective function -----
0064 %% polynomial cost of P and Q
0065 % use totcost only on polynomial cost; in the minimization problem
0066 % formulation, pwl cost is the sum of the y variables.
0067 ipol = find(gencost(:, MODEL) == POLYNOMIAL);   %% poly MW and MVAr costs
0068 xx = [ Pg; Qg ] * baseMVA;
0069 if ~isempty(ipol)
0070   f = sum( totcost(gencost(ipol, :), xx(ipol)) );  %% cost of poly P or Q
0071 else
0072   f = 0;
0073 end
0074 
0075 %% piecewise linear cost of P and Q
0076 if ny > 0
0077   ccost = full(sparse(ones(1,ny), vv.i1.y:vv.iN.y, ones(1,ny), 1, nxyz));
0078   f = f + ccost * x;
0079 else
0080   ccost = zeros(1, nxyz);
0081 end
0082 
0083 %% generalized cost term
0084 if ~isempty(N)
0085     nw = size(N, 1);
0086     r = N * x - rh;                 %% Nx - rhat
0087     iLT = find(r < -kk);            %% below dead zone
0088     iEQ = find(r == 0 & kk == 0);   %% dead zone doesn't exist
0089     iGT = find(r > kk);             %% above dead zone
0090     iND = [iLT; iEQ; iGT];          %% rows that are Not in the Dead region
0091     iL = find(dd == 1);             %% rows using linear function
0092     iQ = find(dd == 2);             %% rows using quadratic function
0093     LL = sparse(iL, iL, 1, nw, nw);
0094     QQ = sparse(iQ, iQ, 1, nw, nw);
0095     kbar = sparse(iND, iND, [   ones(length(iLT), 1);
0096                                 zeros(length(iEQ), 1);
0097                                 -ones(length(iGT), 1)], nw, nw) * kk;
0098     rr = r + kbar;                  %% apply non-dead zone shift
0099     M = sparse(iND, iND, mm(iND), nw, nw);  %% dead zone or scale
0100     diagrr = sparse(1:nw, 1:nw, rr, nw, nw);
0101     
0102     %% linear rows multiplied by rr(i), quadratic rows by rr(i)^2
0103     w = M * (LL + QQ * diagrr) * rr;
0104 
0105     f = f + (w' * H * w) / 2 + Cw' * w;
0106 end
0107 
0108 %%----- evaluate cost gradient -----
0109 if nargout > 1
0110   %% index ranges
0111   iPg = vv.i1.Pg:vv.iN.Pg;
0112   iQg = vv.i1.Qg:vv.iN.Qg;
0113 
0114   %% polynomial cost of P and Q
0115   df_dPgQg = zeros(2*ng, 1);        %% w.r.t p.u. Pg and Qg
0116   df_dPgQg(ipol) = baseMVA * polycost(gencost(ipol, :), xx(ipol), 1);
0117   df = zeros(nxyz, 1);
0118   df(iPg) = df_dPgQg(1:ng);
0119   df(iQg) = df_dPgQg((1:ng) + ng);
0120 
0121   %% piecewise linear cost of P and Q
0122   df = df + ccost';  % As in MINOS, the linear cost row is additive wrt
0123                      % any nonlinear cost.
0124 
0125   %% generalized cost term
0126   if ~isempty(N)
0127     HwC = H * w + Cw;
0128     AA = N' * M * (LL + 2 * QQ * diagrr);
0129     df = df + AA * HwC;
0130     
0131     %% numerical check
0132     if 0    %% 1 to check, 0 to skip check
0133       ddff = zeros(size(df));
0134       step = 1e-7;
0135       tol  = 1e-3;
0136       for k = 1:length(x)
0137         xx = x;
0138         xx(k) = xx(k) + step;
0139         ddff(k) = (opf_costfcn(xx, om) - f) / step;
0140       end
0141       if max(abs(ddff - df)) > tol
0142         idx = find(abs(ddff - df) == max(abs(ddff - df)));
0143         fprintf('\nMismatch in gradient\n');
0144         fprintf('idx             df(num)         df              diff\n');
0145         fprintf('%4d%16g%16g%16g\n', [ 1:length(df); ddff'; df'; abs(ddff - df)' ]);
0146         fprintf('MAX\n');
0147         fprintf('%4d%16g%16g%16g\n', [ idx'; ddff(idx)'; df(idx)'; abs(ddff(idx) - df(idx))' ]);
0148         fprintf('\n');
0149       end
0150     end     %% numeric check
0151   end
0152 
0153   %% ---- evaluate cost Hessian -----
0154   if nargout > 2
0155     pcost = gencost(1:ng, :);
0156     if size(gencost, 1) > ng
0157         qcost = gencost(ng+1:2*ng, :);
0158     else
0159         qcost = [];
0160     end
0161     
0162     %% polynomial generator costs
0163     d2f_dPg2 = sparse(ng, 1);               %% w.r.t. p.u. Pg
0164     d2f_dQg2 = sparse(ng, 1);               %% w.r.t. p.u. Qg
0165     ipolp = find(pcost(:, MODEL) == POLYNOMIAL);
0166     d2f_dPg2(ipolp) = baseMVA^2 * polycost(pcost(ipolp, :), Pg(ipolp)*baseMVA, 2);
0167     if ~isempty(qcost)          %% Qg is not free
0168         ipolq = find(qcost(:, MODEL) == POLYNOMIAL);
0169         d2f_dQg2(ipolq) = baseMVA^2 * polycost(qcost(ipolq, :), Qg(ipolq)*baseMVA, 2);
0170     end
0171     i = [iPg iQg]';
0172     d2f = sparse(i, i, [d2f_dPg2; d2f_dQg2], nxyz, nxyz);
0173 
0174     %% generalized cost
0175     if ~isempty(N)
0176         d2f = d2f + AA * H * AA' + 2 * N' * M * QQ * sparse(1:nw, 1:nw, HwC, nw, nw) * N;
0177     end
0178   end
0179 end

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