Home > matpower6.0 > t > t_mips.m

t_mips

PURPOSE ^

T_MIPS Tests of MIPS NLP solver.

SYNOPSIS ^

function t_mips(quiet)

DESCRIPTION ^

T_MIPS  Tests of MIPS NLP solver.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SUBFUNCTIONS ^

SOURCE CODE ^

0001 function t_mips(quiet)
0002 %T_MIPS  Tests of MIPS NLP solver.
0003 
0004 %   MIPS
0005 %   Copyright (c) 2010-2016, Power Systems Engineering Research Center (PSERC)
0006 %   by Ray Zimmerman, PSERC Cornell
0007 %
0008 %   This file is part of MIPS.
0009 %   Covered by the 3-clause BSD License (see LICENSE file for details).
0010 %   See http://www.pserc.cornell.edu/matpower/ for more info.
0011 
0012 if nargin < 1
0013     quiet = 0;
0014 end
0015 
0016 t_begin(60, quiet);
0017 
0018 t = 'unconstrained banana function : ';
0019 %% from MATLAB Optimization Toolbox's bandem.m
0020 f_fcn = @(x)f2(x);
0021 x0 = [-1.9; 2];
0022 % [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], [], [], [], [], struct('verbose', 2));
0023 [x, f, s, out, lam] = mips(f_fcn, x0);
0024 t_is(s, 1, 13, [t 'success']);
0025 t_is(x, [1; 1], 13, [t 'x']);
0026 t_is(f, 0, 13, [t 'f']);
0027 t_is(out.hist(end).compcond, 0, 6, [t 'compcond']);
0028 t_ok(isempty(lam.mu_l), [t 'lam.mu_l']);
0029 t_ok(isempty(lam.mu_u), [t 'lam.mu_u']);
0030 t_is(lam.lower, zeros(size(x)), 13, [t 'lam.lower']);
0031 t_is(lam.upper, zeros(size(x)), 13, [t 'lam.upper']);
0032 
0033 t = 'unconstrained 3-d quadratic : ';
0034 %% from http://www.akiti.ca/QuadProgEx0Constr.html
0035 f_fcn = @(x)f3(x);
0036 x0 = [0; 0; 0];
0037 % [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], [], [], [], [], struct('verbose', 2));
0038 [x, f, s, out, lam] = mips(f_fcn, x0);
0039 t_is(s, 1, 13, [t 'success']);
0040 t_is(x, [3; 5; 7], 13, [t 'x']);
0041 t_is(f, -244, 13, [t 'f']);
0042 t_is(out.hist(end).compcond, 0, 6, [t 'compcond']);
0043 t_ok(isempty(lam.mu_l), [t 'lam.mu_l']);
0044 t_ok(isempty(lam.mu_u), [t 'lam.mu_u']);
0045 t_is(lam.lower, zeros(size(x)), 13, [t 'lam.lower']);
0046 t_is(lam.upper, zeros(size(x)), 13, [t 'lam.upper']);
0047 
0048 t = 'constrained 4-d QP : ';
0049 %% from http://www.jmu.edu/docs/sasdoc/sashtml/iml/chap8/sect12.htm
0050 f_fcn = @(x)f4(x);
0051 x0 = [1; 0; 0; 1];
0052 A = [   1       1       1       1;
0053         0.17    0.11    0.10    0.18    ];
0054 l = [1; 0.10];
0055 u = [1; Inf];
0056 xmin = zeros(4,1);
0057 % [x, f, s, out, lam] = mips(f_fcn, x0, A, l, u, xmin, [], [], [], struct('verbose', 2));
0058 [x, f, s, out, lam] = mips(f_fcn, x0, A, l, u, xmin);
0059 t_is(s, 1, 13, [t 'success']);
0060 t_is(x, [0; 2.8; 0.2; 0]/3, 6, [t 'x']);
0061 t_is(f, 3.29/3, 6, [t 'f']);
0062 t_is(out.hist(end).compcond, 0, 6, [t 'compcond']);
0063 t_is(lam.mu_l, [6.58;0]/3, 6, [t 'lam.mu_l']);
0064 t_is(lam.mu_u, [0;0], 13, [t 'lam.mu_u']);
0065 t_is(lam.lower, [2.24;0;0;1.7667], 4, [t 'lam.lower']);
0066 t_is(lam.upper, zeros(size(x)), 13, [t 'lam.upper']);
0067 
0068 H = [   1003.1  4.3     6.3     5.9;
0069         4.3     2.2     2.1     3.9;
0070         6.3     2.1     3.5     4.8;
0071         5.9     3.9     4.8     10  ];
0072 c = zeros(4,1);
0073 % %% check with quadprog (for dev testing only)
0074 % [x, f, s, out, lam] = quadprog(H,c,-A(2,:), -0.10, A(1,:), 1, xmin);
0075 % t_is(s, 1, 13, [t 'success']);
0076 % t_is(x, [0; 2.8; 0.2; 0]/3, 6, [t 'x']);
0077 % t_is(f, 3.29/3, 6, [t 'f']);
0078 % t_is(lam.eqlin, -6.58/3, 6, [t 'lam.eqlin']);
0079 % t_is(lam.ineqlin, 0, 13, [t 'lam.ineqlin']);
0080 % t_is(lam.lower, [2.24;0;0;1.7667], 4, [t 'lam.lower']);
0081 % t_is(lam.upper, [0;0;0;0], 13, [t 'lam.upper']);
0082 
0083 t = 'constrained 2-d nonlinear : ';
0084 %% from http://en.wikipedia.org/wiki/Nonlinear_programming#2-dimensional_example
0085 f_fcn = @(x)f5(x);
0086 gh_fcn = @(x)gh5(x);
0087 hess_fcn = @(x, lam, cost_mult)hess5(x, lam, cost_mult);
0088 x0 = [1.1; 0];
0089 xmin = zeros(2, 1);
0090 % xmax = 3 * ones(2, 1);
0091 % [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], xmin, [], gh_fcn, hess_fcn, struct('verbose', 2));
0092 [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], xmin, [], gh_fcn, hess_fcn);
0093 t_is(s, 1, 13, [t 'success']);
0094 t_is(x, [1; 1], 6, [t 'x']);
0095 t_is(f, -2, 6, [t 'f']);
0096 t_is(out.hist(end).compcond, 0, 6, [t 'compcond']);
0097 t_is(lam.ineqnonlin, [0;0.5], 6, [t 'lam.ineqnonlin']);
0098 t_ok(isempty(lam.mu_l), [t 'lam.mu_l']);
0099 t_ok(isempty(lam.mu_u), [t 'lam.mu_u']);
0100 t_is(lam.lower, zeros(size(x)), 13, [t 'lam.lower']);
0101 t_is(lam.upper, zeros(size(x)), 13, [t 'lam.upper']);
0102 % %% check with fmincon (for dev testing only)
0103 % % fmoptions = optimset('Algorithm', 'interior-point');
0104 % % [x, f, s, out, lam] = fmincon(f_fcn, x0, [], [], [], [], xmin, [], gh_fcn, fmoptions);
0105 % [x, f, s, out, lam] = fmincon(f_fcn, x0, [], [], [], [], [], [], gh_fcn);
0106 % t_is(s, 1, 13, [t 'success']);
0107 % t_is(x, [1; 1], 4, [t 'x']);
0108 % t_is(f, -2, 6, [t 'f']);
0109 % t_is(lam.ineqnonlin, [0;0.5], 6, [t 'lam.ineqnonlin']);
0110 
0111 t = 'constrained 3-d nonlinear : ';
0112 %% from http://en.wikipedia.org/wiki/Nonlinear_programming#3-dimensional_example
0113 f_fcn = @(x)f6(x);
0114 gh_fcn = @(x)gh6(x);
0115 hess_fcn = @(x, lam, cost_mult)hess6(x, lam, cost_mult);
0116 x0 = [1; 1; 0];
0117 % [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], [], [], gh_fcn, hess_fcn, struct('verbose', 2, 'comptol', 1e-9));
0118 [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], [], [], gh_fcn, hess_fcn);
0119 t_is(s, 1, 13, [t 'success']);
0120 t_is(x, [1.58113883; 2.23606798; 1.58113883], 6, [t 'x']);
0121 t_is(f, -5*sqrt(2), 6, [t 'f']);
0122 t_is(out.hist(end).compcond, 0, 6, [t 'compcond']);
0123 t_is(lam.ineqnonlin, [0;sqrt(2)/2], 7, [t 'lam.ineqnonlin']);
0124 t_ok(isempty(lam.mu_l), [t 'lam.mu_l']);
0125 t_ok(isempty(lam.mu_u), [t 'lam.mu_u']);
0126 t_is(lam.lower, zeros(size(x)), 13, [t 'lam.lower']);
0127 t_is(lam.upper, zeros(size(x)), 13, [t 'lam.upper']);
0128 % %% check with fmincon (for dev testing only)
0129 % % fmoptions = optimset('Algorithm', 'interior-point');
0130 % % [x, f, s, out, lam] = fmincon(f_fcn, x0, [], [], [], [], xmin, [], gh_fcn, fmoptions);
0131 % [x, f, s, out, lam] = fmincon(f_fcn, x0, [], [], [], [], [], [], gh_fcn);
0132 % t_is(s, 1, 13, [t 'success']);
0133 % t_is(x, [1.58113883; 2.23606798; 1.58113883], 4, [t 'x']);
0134 % t_is(f, -5*sqrt(2), 8, [t 'f']);
0135 % t_is(lam.ineqnonlin, [0;sqrt(2)/2], 8, [t 'lam.ineqnonlin']);
0136 
0137 t = 'constrained 3-d nonlinear (struct) : ';
0138 p = struct('f_fcn', f_fcn, 'x0', x0, 'gh_fcn', gh_fcn, 'hess_fcn', hess_fcn);
0139 [x, f, s, out, lam] = mips(p);
0140 t_is(s, 1, 13, [t 'success']);
0141 t_is(x, [1.58113883; 2.23606798; 1.58113883], 6, [t 'x']);
0142 t_is(f, -5*sqrt(2), 6, [t 'f']);
0143 t_is(out.hist(end).compcond, 0, 6, [t 'compcond']);
0144 t_is(lam.ineqnonlin, [0;sqrt(2)/2], 7, [t 'lam.ineqnonlin']);
0145 t_ok(isempty(lam.mu_l), [t 'lam.mu_l']);
0146 t_ok(isempty(lam.mu_u), [t 'lam.mu_u']);
0147 t_is(lam.lower, zeros(size(x)), 13, [t 'lam.lower']);
0148 t_is(lam.upper, zeros(size(x)), 13, [t 'lam.upper']);
0149 
0150 t = 'constrained 4-d nonlinear : ';
0151 %% Hock & Schittkowski test problem #71
0152 f_fcn = @(x)f7(x);
0153 gh_fcn = @(x)gh7(x);
0154 hess_fcn = @(x, lam, sigma)hess7(x, lam, sigma);
0155 x0 = [1; 5; 5; 1];
0156 xmin = ones(4, 1);
0157 xmax = 5 * xmin;
0158 % [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], xmin, xmax, gh_fcn, hess_fcn, struct('verbose', 2, 'comptol', 1e-9));
0159 [x, f, s, out, lam] = mips(f_fcn, x0, [], [], [], xmin, xmax, gh_fcn, hess_fcn);
0160 t_is(s, 1, 13, [t 'success']);
0161 t_is(x, [1; 4.7429994; 3.8211503; 1.3794082], 6, [t 'x']);
0162 t_is(f, 17.0140173, 6, [t 'f']);
0163 t_is(lam.eqnonlin, 0.1614686, 5, [t 'lam.eqnonlin']);
0164 t_is(lam.ineqnonlin, 0.55229366, 5, [t 'lam.ineqnonlin']);
0165 t_ok(isempty(lam.mu_l), [t 'lam.mu_l']);
0166 t_ok(isempty(lam.mu_u), [t 'lam.mu_u']);
0167 t_is(lam.lower, [1.08787121024; 0; 0; 0], 5, [t 'lam.lower']);
0168 t_is(lam.upper, zeros(size(x)), 7, [t 'lam.upper']);
0169 
0170 t_end;
0171 
0172 
0173 % %%-----  eg99 : linearly constrained fmincon example, mips can't solve  -----
0174 % function [f, df, d2f] = eg99(x)
0175 % f = -x(1)*x(2)*x(3);
0176 % df = -[ x(2)*x(3);
0177 %         x(1)*x(3);
0178 %         x(1)*x(2)   ];
0179 % d2f = -[    0       x(3)    x(2);
0180 %             x(3)    0       x(1);
0181 %             x(2)    x(1)    0   ];
0182 % end
0183 %
0184 % x0 = [10;10;10];
0185 % A = [1 2 2];
0186 % l = 0;
0187 % u = 72;
0188 % fmoptions = optimset('Display', 'testing');
0189 % fmoptions = optimset(fmoptions, 'Algorithm', 'interior-point');
0190 % [x, f, s, out, lam] = fmincon(f_fcn, x0, [-A; A], [-l; u], [], [], [], [], [], fmoptions);
0191 % t_is(x, [24; 12; 12], 13, t);
0192 % t_is(f, -3456, 13, t);
0193 
0194 
0195 %% unconstrained banana function
0196 %% from MATLAB Optimization Toolbox's bandem.m
0197 function [f, df, d2f] = f2(x)
0198     a = 100;
0199     f = a*(x(2)-x(1)^2)^2+(1-x(1))^2;
0200     df = [  4*a*(x(1)^3 - x(1)*x(2)) + 2*x(1)-2;
0201             2*a*(x(2) - x(1)^2)                     ];
0202     d2f = 4*a*[ 3*x(1)^2 - x(2) + 1/(2*a),  -x(1);
0203                 -x(1)                       1/2       ];
0204 
0205 
0206 %% unconstrained 3-d quadratic
0207 %% from http://www.akiti.ca/QuadProgEx0Constr.html
0208 function [f, df, d2f] = f3(x)
0209     H = [5 -2 -1; -2 4 3; -1 3 5];
0210     c = [2; -35; -47];
0211     f = 1/2 * x'*H*x + c'*x + 5;
0212     df = H*x + c;
0213     d2f = H;
0214 
0215 
0216 %% constrained 4-d QP
0217 %% from http://www.jmu.edu/docs/sasdoc/sashtml/iml/chap8/sect12.htm
0218 function [f, df, d2f] = f4(x)
0219     H = [   1003.1  4.3     6.3     5.9;
0220             4.3     2.2     2.1     3.9;
0221             6.3     2.1     3.5     4.8;
0222             5.9     3.9     4.8     10  ];
0223     c = zeros(4,1);
0224     f = 1/2 * x'*H*x + c'*x;
0225     df = H*x + c;
0226     d2f = H;
0227 
0228 
0229 %% constrained 2-d nonlinear
0230 %% from http://en.wikipedia.org/wiki/Nonlinear_programming#2-dimensional_example
0231 function [f, df, d2f] = f5(x)
0232     c = -[1; 1];
0233     f = c'*x;
0234     df = c;
0235     d2f = zeros(2,2);
0236 
0237 function [h, g, dh, dg] = gh5(x)
0238     h = [ -1 -1; 1 1] * x.^2 + [1; -2];
0239     dh = 2 * [-x(1) x(1); -x(2) x(2)];
0240     g = []; dg = [];
0241 
0242 function Lxx = hess5(x, lam, cost_mult)
0243     mu = lam.ineqnonlin;
0244     Lxx = 2*[-1 1]*mu*eye(2);
0245 
0246 
0247 %% constrained 3-d nonlinear
0248 %% from http://en.wikipedia.org/wiki/Nonlinear_programming#3-dimensional_example
0249 function [f, df, d2f] = f6(x)
0250     f = -x(1)*x(2) - x(2)*x(3);
0251     df = -[x(2); x(1)+x(3); x(2)];
0252     d2f = -[0 1 0; 1 0 1; 0 1 0];
0253 
0254 function [h, g, dh, dg] = gh6(x)
0255     h = [ 1 -1 1; 1 1 1] * x.^2 + [-2; -10];
0256     dh = 2 * [x(1) x(1); -x(2) x(2); x(3) x(3)];
0257     g = []; dg = [];
0258 
0259 function Lxx = hess6(x, lam, cost_mult)
0260     if nargin < 3, cost_mult = 1; end
0261     mu = lam.ineqnonlin;
0262     Lxx = cost_mult * [0 -1 0; -1 0 -1; 0 -1 0] + ...
0263             [2*[1 1]*mu 0 0; 0 2*[-1 1]*mu 0; 0 0 2*[1 1]*mu];
0264 
0265 
0266 %% constrained 4-d nonlinear
0267 %% Hock & Schittkowski test problem #71
0268 function [f, df, d2f] = f7(x)
0269     f = x(1)*x(4)*sum(x(1:3)) + x(3);
0270     df = [ x(1)*x(4) + x(4)*sum(x(1:3));
0271            x(1)*x(4);
0272            x(1)*x(4) + 1;
0273            x(1)*sum(x(1:3)) ];
0274     d2f = sparse([ 2*x(4)        x(4)   x(4)  2*x(1)+x(2)+x(3);
0275               x(4)               0      0     x(1);
0276               x(4)               0      0     x(1);
0277               2*x(1)+x(2)+x(3)  x(1)  x(1)    0
0278         ]);
0279 
0280 function [h, g, dh, dg] = gh7(x)
0281     g = sum(x.^2) - 40;
0282     h = -prod(x) + 25;
0283     dg = 2*x;
0284     dh = -prod(x)./x;
0285 
0286 function Lxx = hess7(x, lam, sigma)
0287     if nargin < 3, sigma = 1; end
0288     lambda = lam.eqnonlin;
0289     mu     = lam.ineqnonlin;
0290     [f, df, d2f] = f7(x);
0291     Lxx = sigma * d2f + lambda*2*speye(4) - ...
0292        mu*sparse([      0     x(3)*x(4) x(2)*x(4) x(2)*x(3);
0293                     x(3)*x(4)     0     x(1)*x(4) x(1)*x(3);
0294                     x(2)*x(4) x(1)*x(4)     0     x(1)*x(2);
0295                     x(2)*x(3) x(1)*x(3) x(1)*x(2)     0  ]);

Generated on Fri 16-Dec-2016 12:45:37 by m2html © 2005