Home > matpower4.0 > fmincopf_solver.m

fmincopf_solver

PURPOSE ^

FMINCOPF_SOLVER Solves AC optimal power flow using FMINCON.

SYNOPSIS ^

function [results, success, raw] = fmincopf_solver(om, mpopt)

DESCRIPTION ^

FMINCOPF_SOLVER  Solves AC optimal power flow using FMINCON.

   [RESULTS, SUCCESS, RAW] = FMINCOPF_SOLVER(OM, MPOPT)

   Inputs are an OPF model object and a MATPOWER options vector.

   Outputs are a RESULTS struct, SUCCESS flag and RAW output struct.

   RESULTS is a MATPOWER case struct (mpc) with the usual baseMVA, bus
   branch, gen, gencost fields, along with the following additional
   fields:
       .order      see 'help ext2int' for details of this field
       .x          final value of optimization variables (internal order)
       .f          final objective function value
       .mu         shadow prices on ...
           .var
               .l  lower bounds on variables
               .u  upper bounds on variables
           .nln
               .l  lower bounds on nonlinear constraints
               .u  upper bounds on nonlinear constraints
           .lin
               .l  lower bounds on linear constraints
               .u  upper bounds on linear constraints

   SUCCESS     1 if solver converged successfully, 0 otherwise

   RAW         raw output in form returned by MINOS
       .xr     final value of optimization variables
       .pimul  constraint multipliers
       .info   solver specific termination code
       .output solver specific output information

   See also OPF, FMINCON.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [results, success, raw] = fmincopf_solver(om, mpopt)
0002 %FMINCOPF_SOLVER  Solves AC optimal power flow using FMINCON.
0003 %
0004 %   [RESULTS, SUCCESS, RAW] = FMINCOPF_SOLVER(OM, MPOPT)
0005 %
0006 %   Inputs are an OPF model object and a MATPOWER options vector.
0007 %
0008 %   Outputs are a RESULTS struct, SUCCESS flag and RAW output struct.
0009 %
0010 %   RESULTS is a MATPOWER case struct (mpc) with the usual baseMVA, bus
0011 %   branch, gen, gencost fields, along with the following additional
0012 %   fields:
0013 %       .order      see 'help ext2int' for details of this field
0014 %       .x          final value of optimization variables (internal order)
0015 %       .f          final objective function value
0016 %       .mu         shadow prices on ...
0017 %           .var
0018 %               .l  lower bounds on variables
0019 %               .u  upper bounds on variables
0020 %           .nln
0021 %               .l  lower bounds on nonlinear constraints
0022 %               .u  upper bounds on nonlinear constraints
0023 %           .lin
0024 %               .l  lower bounds on linear constraints
0025 %               .u  upper bounds on linear constraints
0026 %
0027 %   SUCCESS     1 if solver converged successfully, 0 otherwise
0028 %
0029 %   RAW         raw output in form returned by MINOS
0030 %       .xr     final value of optimization variables
0031 %       .pimul  constraint multipliers
0032 %       .info   solver specific termination code
0033 %       .output solver specific output information
0034 %
0035 %   See also OPF, FMINCON.
0036 
0037 %   MATPOWER
0038 %   $Id: fmincopf_solver.m,v 1.23 2010/06/09 14:56:58 ray Exp $
0039 %   by Ray Zimmerman, PSERC Cornell
0040 %   and Carlos E. Murillo-Sanchez, PSERC Cornell & Universidad Autonoma de Manizales
0041 %   Copyright (c) 2000-2010 by Power System Engineering Research Center (PSERC)
0042 %
0043 %   This file is part of MATPOWER.
0044 %   See http://www.pserc.cornell.edu/matpower/ for more info.
0045 %
0046 %   MATPOWER is free software: you can redistribute it and/or modify
0047 %   it under the terms of the GNU General Public License as published
0048 %   by the Free Software Foundation, either version 3 of the License,
0049 %   or (at your option) any later version.
0050 %
0051 %   MATPOWER is distributed in the hope that it will be useful,
0052 %   but WITHOUT ANY WARRANTY; without even the implied warranty of
0053 %   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
0054 %   GNU General Public License for more details.
0055 %
0056 %   You should have received a copy of the GNU General Public License
0057 %   along with MATPOWER. If not, see <http://www.gnu.org/licenses/>.
0058 %
0059 %   Additional permission under GNU GPL version 3 section 7
0060 %
0061 %   If you modify MATPOWER, or any covered work, to interface with
0062 %   other modules (such as MATLAB code and MEX-files) available in a
0063 %   MATLAB(R) or comparable environment containing parts covered
0064 %   under other licensing terms, the licensors of MATPOWER grant
0065 %   you additional permission to convey the resulting work.
0066 
0067 %%----- initialization -----
0068 %% define named indices into data matrices
0069 [PQ, PV, REF, NONE, BUS_I, BUS_TYPE, PD, QD, GS, BS, BUS_AREA, VM, ...
0070     VA, BASE_KV, ZONE, VMAX, VMIN, LAM_P, LAM_Q, MU_VMAX, MU_VMIN] = idx_bus;
0071 [GEN_BUS, PG, QG, QMAX, QMIN, VG, MBASE, GEN_STATUS, PMAX, PMIN, ...
0072     MU_PMAX, MU_PMIN, MU_QMAX, MU_QMIN, PC1, PC2, QC1MIN, QC1MAX, ...
0073     QC2MIN, QC2MAX, RAMP_AGC, RAMP_10, RAMP_30, RAMP_Q, APF] = idx_gen;
0074 [F_BUS, T_BUS, BR_R, BR_X, BR_B, RATE_A, RATE_B, RATE_C, ...
0075     TAP, SHIFT, BR_STATUS, PF, QF, PT, QT, MU_SF, MU_ST, ...
0076     ANGMIN, ANGMAX, MU_ANGMIN, MU_ANGMAX] = idx_brch;
0077 
0078 %% options
0079 verbose = mpopt(31);    %% VERBOSE
0080 
0081 %% unpack data
0082 mpc = get_mpc(om);
0083 [baseMVA, bus, gen, branch] = ...
0084     deal(mpc.baseMVA, mpc.bus, mpc.gen, mpc.branch);
0085 [vv, ll, nn] = get_idx(om);
0086 
0087 %% problem dimensions
0088 nb = size(bus, 1);          %% number of buses
0089 nl = size(branch, 1);       %% number of branches
0090 ny = getN(om, 'var', 'y');  %% number of piece-wise linear costs
0091 
0092 %% bounds on optimization vars
0093 [x0, LB, UB] = getv(om);
0094 
0095 %% linear constraints
0096 [A, l, u] = linear_constraints(om);
0097 
0098 %% so, can we do anything good about lambda initialization?
0099 if all(bus(:, LAM_P) == 0)
0100   bus(:, LAM_P) = (10)*ones(nb, 1);
0101 end
0102 
0103 %% split l <= A*x <= u into less than, equal to, greater than, and
0104 %% doubly-bounded sets
0105 ieq = find( abs(u-l) <= eps );          %% equality
0106 igt = find( u >=  1e10 & l > -1e10 );   %% greater than, unbounded above
0107 ilt = find( l <= -1e10 & u <  1e10 );   %% less than, unbounded below
0108 ibx = find( (abs(u-l) > eps) & (u < 1e10) & (l > -1e10) );
0109 Af  = [ A(ilt, :); -A(igt, :); A(ibx, :); -A(ibx, :) ];
0110 bf  = [ u(ilt);   -l(igt);     u(ibx);    -l(ibx)];
0111 Afeq = A(ieq, :);
0112 bfeq = u(ieq);
0113 
0114 %% build admittance matrices
0115 [Ybus, Yf, Yt] = makeYbus(baseMVA, bus, branch);
0116 
0117 %% tolerances
0118 if mpopt(19) == 0           %% CONSTR_MAX_IT
0119   mpopt(19) = 150 + 2*nb;
0120 end
0121 
0122 %% find branches with flow limits
0123 il = find(branch(:, RATE_A) ~= 0 & branch(:, RATE_A) < 1e10);
0124 nl2 = length(il);           %% number of constrained lines
0125 
0126 %% basic optimset options needed for fmincon
0127 fmoptions = optimset('GradObj', 'on', 'GradConstr', 'on', ...
0128             'MaxIter', mpopt(19), 'TolCon', mpopt(16), ...
0129             'TolX', mpopt(17), 'TolFun', mpopt(18) );
0130 fmoptions.MaxFunEvals = 4 * fmoptions.MaxIter;
0131 if verbose == 0,
0132   fmoptions.Display = 'off';
0133 elseif verbose == 1
0134   fmoptions.Display = 'iter';
0135 else
0136   fmoptions.Display = 'testing';
0137 end
0138 
0139 %% select algorithm
0140 otver = ver('optim');
0141 if str2double(otver.Version(1)) < 4
0142   fmoptions = optimset(fmoptions, 'LargeScale', 'off');
0143   Af = full(Af);
0144   Afeq = full(Afeq);
0145 else
0146   if mpopt(55) == 1           %% active-set
0147     fmoptions = optimset(fmoptions, 'Algorithm', 'active-set');
0148     Af = full(Af);
0149     Afeq = full(Afeq);
0150   elseif mpopt(55) == 2       %% interior-point, w/ default 'bfgs' Hessian approx
0151     fmoptions = optimset(fmoptions, 'Algorithm', 'interior-point');
0152   elseif mpopt(55) == 3       %% interior-point, w/ 'lbfgs' Hessian approx
0153     fmoptions = optimset(fmoptions, 'Algorithm', 'interior-point', 'Hessian','lbfgs');
0154   elseif mpopt(55) == 4       %% interior-point, w/ exact user-supplied Hessian
0155     fmc_hessian = @(x, lambda)opf_hessfcn(x, lambda, 1, om, Ybus, Yf(il,:), Yt(il,:), mpopt, il);
0156     fmoptions = optimset(fmoptions, 'Algorithm', 'interior-point', ...
0157         'Hessian', 'user-supplied', 'HessFcn', fmc_hessian);
0158   elseif mpopt(55) == 5       %% interior-point, w/ finite-diff Hessian
0159     fmoptions = optimset(fmoptions, 'Algorithm', 'interior-point', 'Hessian','fin-diff-grads', 'SubProblem', 'cg');
0160   else
0161     error('fmincopf_solver: unknown algorithm specified in FMC_ALG option');
0162   end
0163 end
0164 % fmoptions = optimset(fmoptions, 'DerivativeCheck', 'on', 'FinDiffType', 'central', 'FunValCheck', 'on');
0165 % fmoptions = optimset(fmoptions, 'Diagnostics', 'on');
0166 
0167 %% try to select an interior initial point for interior point solver
0168 if str2double(otver.Version(1)) >= 4 && strcmp(optimget(fmoptions, 'Algorithm'), 'interior-point')
0169   x0 = zeros(getN(om, 'var'), 1);
0170   x0(vv.i1.Va:vv.iN.Va) = 0;
0171   x0(vv.i1.Vm:vv.iN.Vm)   = 1;
0172   x0(vv.i1.Pg:vv.iN.Pg) = (gen(:, PMIN) + gen(:, PMAX)) / 2 / baseMVA;
0173   x0(vv.i1.Qg:vv.iN.Qg) = (gen(:, QMIN) + gen(:, QMAX)) / 2 / baseMVA;
0174 end
0175 
0176 %%-----  run opf  -----
0177 f_fcn = @(x)opf_costfcn(x, om);
0178 gh_fcn = @(x)opf_consfcn(x, om, Ybus, Yf(il,:), Yt(il,:), mpopt, il);
0179 [x, f, info, Output, Lambda] = ...
0180   fmincon(f_fcn, x0, Af, bf, Afeq, bfeq, LB, UB, gh_fcn, fmoptions);
0181 success = (info > 0);
0182 
0183 %% update solution data
0184 Va = x(vv.i1.Va:vv.iN.Va);
0185 Vm = x(vv.i1.Vm:vv.iN.Vm);
0186 Pg = x(vv.i1.Pg:vv.iN.Pg);
0187 Qg = x(vv.i1.Qg:vv.iN.Qg);
0188 V = Vm .* exp(1j*Va);
0189 
0190 %%-----  calculate return values  -----
0191 %% update voltages & generator outputs
0192 bus(:, VA) = Va * 180/pi;
0193 bus(:, VM) = Vm;
0194 gen(:, PG) = Pg * baseMVA;
0195 gen(:, QG) = Qg * baseMVA;
0196 gen(:, VG) = Vm(gen(:, GEN_BUS));
0197 
0198 %% compute branch flows
0199 Sf = V(branch(:, F_BUS)) .* conj(Yf * V);  %% cplx pwr at "from" bus, p.u.
0200 St = V(branch(:, T_BUS)) .* conj(Yt * V);  %% cplx pwr at "to" bus, p.u.
0201 branch(:, PF) = real(Sf) * baseMVA;
0202 branch(:, QF) = imag(Sf) * baseMVA;
0203 branch(:, PT) = real(St) * baseMVA;
0204 branch(:, QT) = imag(St) * baseMVA;
0205 
0206 %% line constraint is actually on square of limit
0207 %% so we must fix multipliers
0208 muSf = zeros(nl, 1);
0209 muSt = zeros(nl, 1);
0210 if ~isempty(il)
0211     muSf(il) = 2 * Lambda.ineqnonlin(1:nl2)       .* branch(il, RATE_A) / baseMVA;
0212     muSt(il) = 2 * Lambda.ineqnonlin((1:nl2)+nl2) .* branch(il, RATE_A) / baseMVA;
0213 end
0214 
0215 %% update Lagrange multipliers
0216 bus(:, MU_VMAX)  = Lambda.upper(vv.i1.Vm:vv.iN.Vm);
0217 bus(:, MU_VMIN)  = Lambda.lower(vv.i1.Vm:vv.iN.Vm);
0218 gen(:, MU_PMAX)  = Lambda.upper(vv.i1.Pg:vv.iN.Pg) / baseMVA;
0219 gen(:, MU_PMIN)  = Lambda.lower(vv.i1.Pg:vv.iN.Pg) / baseMVA;
0220 gen(:, MU_QMAX)  = Lambda.upper(vv.i1.Qg:vv.iN.Qg) / baseMVA;
0221 gen(:, MU_QMIN)  = Lambda.lower(vv.i1.Qg:vv.iN.Qg) / baseMVA;
0222 bus(:, LAM_P)    = Lambda.eqnonlin(nn.i1.Pmis:nn.iN.Pmis) / baseMVA;
0223 bus(:, LAM_Q)    = Lambda.eqnonlin(nn.i1.Qmis:nn.iN.Qmis) / baseMVA;
0224 branch(:, MU_SF) = muSf / baseMVA;
0225 branch(:, MU_ST) = muSt / baseMVA;
0226 
0227 %% package up results
0228 nlnN = getN(om, 'nln');
0229 nlt = length(ilt);
0230 ngt = length(igt);
0231 nbx = length(ibx);
0232 
0233 %% extract multipliers for nonlinear constraints
0234 kl = find(Lambda.eqnonlin < 0);
0235 ku = find(Lambda.eqnonlin > 0);
0236 nl_mu_l = zeros(nlnN, 1);
0237 nl_mu_u = [zeros(2*nb, 1); muSf; muSt];
0238 nl_mu_l(kl) = -Lambda.eqnonlin(kl);
0239 nl_mu_u(ku) =  Lambda.eqnonlin(ku);
0240 
0241 %% extract multipliers for linear constraints
0242 kl = find(Lambda.eqlin < 0);
0243 ku = find(Lambda.eqlin > 0);
0244 
0245 mu_l = zeros(size(u));
0246 mu_l(ieq(kl)) = -Lambda.eqlin(kl);
0247 mu_l(igt) = Lambda.ineqlin(nlt+(1:ngt));
0248 mu_l(ibx) = Lambda.ineqlin(nlt+ngt+nbx+(1:nbx));
0249 
0250 mu_u = zeros(size(u));
0251 mu_u(ieq(ku)) = Lambda.eqlin(ku);
0252 mu_u(ilt) = Lambda.ineqlin(1:nlt);
0253 mu_u(ibx) = Lambda.ineqlin(nlt+ngt+(1:nbx));
0254 
0255 mu = struct( ...
0256   'var', struct('l', Lambda.lower, 'u', Lambda.upper), ...
0257   'nln', struct('l', nl_mu_l, 'u', nl_mu_u), ...
0258   'lin', struct('l', mu_l, 'u', mu_u) );
0259 
0260 results = mpc;
0261 [results.bus, results.branch, results.gen, ...
0262     results.om, results.x, results.mu, results.f] = ...
0263         deal(bus, branch, gen, om, x, mu, f);
0264 
0265 pimul = [ ...
0266   results.mu.nln.l - results.mu.nln.u;
0267   results.mu.lin.l - results.mu.lin.u;
0268   -ones(ny>0, 1);
0269   results.mu.var.l - results.mu.var.u;
0270 ];
0271 raw = struct('xr', x, 'pimul', pimul, 'info', info, 'output', Output);

Generated on Mon 26-Jan-2015 14:56:45 by m2html © 2005