annotate pyCSalgos/BP/l1qc.py @ 68:cab8a215f9a1 tip

Minor
author Nic Cleju <nikcleju@gmail.com>
date Tue, 09 Jul 2013 14:50:09 +0300
parents afcfd4d1d548
children
rev   line source
nikcleju@2 1
nikcleju@2 2 import numpy as np
nikcleju@2 3 import scipy.linalg
nikcleju@2 4 import math
nikcleju@2 5
nikcleju@37 6 class l1qcInputValueError(Exception):
nikcleju@37 7 pass
nikcleju@2 8
nikcleju@2 9 #function [x, res, iter] = cgsolve(A, b, tol, maxiter, verbose)
nikcleju@2 10 def cgsolve(A, b, tol, maxiter, verbose=1):
nikcleju@2 11 # Solve a symmetric positive definite system Ax = b via conjugate gradients.
nikcleju@2 12 #
nikcleju@2 13 # Usage: [x, res, iter] = cgsolve(A, b, tol, maxiter, verbose)
nikcleju@2 14 #
nikcleju@2 15 # A - Either an NxN matrix, or a function handle.
nikcleju@2 16 #
nikcleju@2 17 # b - N vector
nikcleju@2 18 #
nikcleju@2 19 # tol - Desired precision. Algorithm terminates when
nikcleju@2 20 # norm(Ax-b)/norm(b) < tol .
nikcleju@2 21 #
nikcleju@2 22 # maxiter - Maximum number of iterations.
nikcleju@2 23 #
nikcleju@2 24 # verbose - If 0, do not print out progress messages.
nikcleju@2 25 # If and integer greater than 0, print out progress every 'verbose' iters.
nikcleju@2 26 #
nikcleju@2 27 # Written by: Justin Romberg, Caltech
nikcleju@2 28 # Email: jrom@acm.caltech.edu
nikcleju@2 29 # Created: October 2005
nikcleju@2 30 #
nikcleju@2 31
nikcleju@2 32 #---------------------
nikcleju@2 33 # Original Matab code:
nikcleju@2 34 #
nikcleju@2 35 #if (nargin < 5), verbose = 1; end
nikcleju@2 36 #
nikcleju@2 37 #implicit = isa(A,'function_handle');
nikcleju@2 38 #
nikcleju@2 39 #x = zeros(length(b),1);
nikcleju@2 40 #r = b;
nikcleju@2 41 #d = r;
nikcleju@2 42 #delta = r'*r;
nikcleju@2 43 #delta0 = b'*b;
nikcleju@2 44 #numiter = 0;
nikcleju@2 45 #bestx = x;
nikcleju@2 46 #bestres = sqrt(delta/delta0);
nikcleju@2 47 #while ((numiter < maxiter) && (delta > tol^2*delta0))
nikcleju@2 48 #
nikcleju@2 49 # # q = A*d
nikcleju@2 50 # if (implicit), q = A(d); else q = A*d; end
nikcleju@2 51 #
nikcleju@2 52 # alpha = delta/(d'*q);
nikcleju@2 53 # x = x + alpha*d;
nikcleju@2 54 #
nikcleju@2 55 # if (mod(numiter+1,50) == 0)
nikcleju@2 56 # # r = b - Aux*x
nikcleju@2 57 # if (implicit), r = b - A(x); else r = b - A*x; end
nikcleju@2 58 # else
nikcleju@2 59 # r = r - alpha*q;
nikcleju@2 60 # end
nikcleju@2 61 #
nikcleju@2 62 # deltaold = delta;
nikcleju@2 63 # delta = r'*r;
nikcleju@2 64 # beta = delta/deltaold;
nikcleju@2 65 # d = r + beta*d;
nikcleju@2 66 # numiter = numiter + 1;
nikcleju@2 67 # if (sqrt(delta/delta0) < bestres)
nikcleju@2 68 # bestx = x;
nikcleju@2 69 # bestres = sqrt(delta/delta0);
nikcleju@2 70 # end
nikcleju@2 71 #
nikcleju@2 72 # if ((verbose) && (mod(numiter,verbose)==0))
nikcleju@2 73 # disp(sprintf('cg: Iter = #d, Best residual = #8.3e, Current residual = #8.3e', ...
nikcleju@2 74 # numiter, bestres, sqrt(delta/delta0)));
nikcleju@2 75 # end
nikcleju@2 76 #
nikcleju@2 77 #end
nikcleju@2 78 #
nikcleju@2 79 #if (verbose)
nikcleju@2 80 # disp(sprintf('cg: Iterations = #d, best residual = #14.8e', numiter, bestres));
nikcleju@2 81 #end
nikcleju@2 82 #x = bestx;
nikcleju@2 83 #res = bestres;
nikcleju@2 84 #iter = numiter;
nikcleju@2 85
nikcleju@2 86 # End of original Matab code
nikcleju@2 87 #----------------------------
nikcleju@2 88
nikcleju@2 89 #if (nargin < 5), verbose = 1; end
nikcleju@2 90 # Optional argument
nikcleju@2 91
nikcleju@2 92 #implicit = isa(A,'function_handle');
nikcleju@2 93 if hasattr(A, '__call__'):
nikcleju@2 94 implicit = True
nikcleju@2 95 else:
nikcleju@2 96 implicit = False
nikcleju@2 97
nikcleju@2 98 x = np.zeros(b.size)
nikcleju@2 99 r = b.copy()
nikcleju@2 100 d = r.copy()
nikcleju@2 101 delta = np.vdot(r,r)
nikcleju@2 102 delta0 = np.vdot(b,b)
nikcleju@2 103 numiter = 0
nikcleju@2 104 bestx = x.copy()
nikcleju@2 105 bestres = math.sqrt(delta/delta0)
nikcleju@2 106 while (numiter < maxiter) and (delta > tol**2*delta0):
nikcleju@2 107
nikcleju@2 108 # q = A*d
nikcleju@2 109 #if (implicit), q = A(d); else q = A*d; end
nikcleju@2 110 if implicit:
nikcleju@2 111 q = A(d)
nikcleju@2 112 else:
nikcleju@2 113 q = np.dot(A,d)
nikcleju@2 114
nikcleju@2 115 alpha = delta/np.vdot(d,q)
nikcleju@2 116 x = x + alpha*d
nikcleju@2 117
nikcleju@2 118 if divmod(numiter+1,50)[1] == 0:
nikcleju@2 119 # r = b - Aux*x
nikcleju@2 120 #if (implicit), r = b - A(x); else r = b - A*x; end
nikcleju@2 121 if implicit:
nikcleju@2 122 r = b - A(x)
nikcleju@2 123 else:
nikcleju@2 124 r = b - np.dot(A,x)
nikcleju@2 125 else:
nikcleju@2 126 r = r - alpha*q
nikcleju@2 127 #end
nikcleju@2 128
nikcleju@2 129 deltaold = delta;
nikcleju@2 130 delta = np.vdot(r,r)
nikcleju@2 131 beta = delta/deltaold;
nikcleju@2 132 d = r + beta*d
nikcleju@2 133 numiter = numiter + 1
nikcleju@2 134 if (math.sqrt(delta/delta0) < bestres):
nikcleju@2 135 bestx = x.copy()
nikcleju@2 136 bestres = math.sqrt(delta/delta0)
nikcleju@2 137 #end
nikcleju@2 138
nikcleju@2 139 if ((verbose) and (divmod(numiter,verbose)[1]==0)):
nikcleju@2 140 #disp(sprintf('cg: Iter = #d, Best residual = #8.3e, Current residual = #8.3e', ...
nikcleju@2 141 # numiter, bestres, sqrt(delta/delta0)));
nikcleju@2 142 print 'cg: Iter = ',numiter,', Best residual = ',bestres,', Current residual = ',math.sqrt(delta/delta0)
nikcleju@2 143 #end
nikcleju@2 144
nikcleju@2 145 #end
nikcleju@2 146
nikcleju@2 147 if (verbose):
nikcleju@2 148 #disp(sprintf('cg: Iterations = #d, best residual = #14.8e', numiter, bestres));
nikcleju@2 149 print 'cg: Iterations = ',numiter,', best residual = ',bestres
nikcleju@2 150 #end
nikcleju@2 151 x = bestx.copy()
nikcleju@2 152 res = bestres
nikcleju@2 153 iter = numiter
nikcleju@2 154
nikcleju@2 155 return x,res,iter
nikcleju@2 156
nikcleju@2 157
nikcleju@2 158
nikcleju@2 159 #function [xp, up, niter] = l1qc_newton(x0, u0, A, At, b, epsilon, tau, newtontol, newtonmaxiter, cgtol, cgmaxiter)
nikcleju@37 160 def l1qc_newton(x0, u0, A, At, b, epsilon, tau, newtontol, newtonmaxiter, cgtol, cgmaxiter, verbose=False):
nikcleju@2 161 # Newton algorithm for log-barrier subproblems for l1 minimization
nikcleju@2 162 # with quadratic constraints.
nikcleju@2 163 #
nikcleju@2 164 # Usage:
nikcleju@2 165 # [xp,up,niter] = l1qc_newton(x0, u0, A, At, b, epsilon, tau,
nikcleju@2 166 # newtontol, newtonmaxiter, cgtol, cgmaxiter)
nikcleju@2 167 #
nikcleju@2 168 # x0,u0 - starting points
nikcleju@2 169 #
nikcleju@2 170 # A - Either a handle to a function that takes a N vector and returns a K
nikcleju@2 171 # vector , or a KxN matrix. If A is a function handle, the algorithm
nikcleju@2 172 # operates in "largescale" mode, solving the Newton systems via the
nikcleju@2 173 # Conjugate Gradients algorithm.
nikcleju@2 174 #
nikcleju@2 175 # At - Handle to a function that takes a K vector and returns an N vector.
nikcleju@2 176 # If A is a KxN matrix, At is ignored.
nikcleju@2 177 #
nikcleju@2 178 # b - Kx1 vector of observations.
nikcleju@2 179 #
nikcleju@2 180 # epsilon - scalar, constraint relaxation parameter
nikcleju@2 181 #
nikcleju@2 182 # tau - Log barrier parameter.
nikcleju@2 183 #
nikcleju@2 184 # newtontol - Terminate when the Newton decrement is <= newtontol.
nikcleju@2 185 # Default = 1e-3.
nikcleju@2 186 #
nikcleju@2 187 # newtonmaxiter - Maximum number of iterations.
nikcleju@2 188 # Default = 50.
nikcleju@2 189 #
nikcleju@2 190 # cgtol - Tolerance for Conjugate Gradients; ignored if A is a matrix.
nikcleju@2 191 # Default = 1e-8.
nikcleju@2 192 #
nikcleju@2 193 # cgmaxiter - Maximum number of iterations for Conjugate Gradients; ignored
nikcleju@2 194 # if A is a matrix.
nikcleju@2 195 # Default = 200.
nikcleju@2 196 #
nikcleju@2 197 # Written by: Justin Romberg, Caltech
nikcleju@2 198 # Email: jrom@acm.caltech.edu
nikcleju@2 199 # Created: October 2005
nikcleju@2 200 #
nikcleju@2 201
nikcleju@2 202 #---------------------
nikcleju@2 203 # Original Matab code:
nikcleju@2 204
nikcleju@3 205 ## check if the mix A is implicit or explicit
nikcleju@2 206 #largescale = isa(A,'function_handle');
nikcleju@2 207 #
nikcleju@2 208 ## line search parameters
nikcleju@2 209 #alpha = 0.01;
nikcleju@2 210 #beta = 0.5;
nikcleju@2 211 #
nikcleju@2 212 #if (~largescale), AtA = A'*A; end
nikcleju@2 213 #
nikcleju@2 214 ## initial point
nikcleju@2 215 #x = x0;
nikcleju@2 216 #u = u0;
nikcleju@2 217 #if (largescale), r = A(x) - b; else r = A*x - b; end
nikcleju@2 218 #fu1 = x - u;
nikcleju@2 219 #fu2 = -x - u;
nikcleju@2 220 #fe = 1/2*(r'*r - epsilon^2);
nikcleju@2 221 #f = sum(u) - (1/tau)*(sum(log(-fu1)) + sum(log(-fu2)) + log(-fe));
nikcleju@2 222 #
nikcleju@2 223 #niter = 0;
nikcleju@2 224 #done = 0;
nikcleju@2 225 #while (~done)
nikcleju@2 226 #
nikcleju@2 227 # if (largescale), atr = At(r); else atr = A'*r; end
nikcleju@2 228 #
nikcleju@2 229 # ntgz = 1./fu1 - 1./fu2 + 1/fe*atr;
nikcleju@2 230 # ntgu = -tau - 1./fu1 - 1./fu2;
nikcleju@2 231 # gradf = -(1/tau)*[ntgz; ntgu];
nikcleju@2 232 #
nikcleju@2 233 # sig11 = 1./fu1.^2 + 1./fu2.^2;
nikcleju@2 234 # sig12 = -1./fu1.^2 + 1./fu2.^2;
nikcleju@2 235 # sigx = sig11 - sig12.^2./sig11;
nikcleju@2 236 #
nikcleju@2 237 # w1p = ntgz - sig12./sig11.*ntgu;
nikcleju@2 238 # if (largescale)
nikcleju@2 239 # h11pfun = @(z) sigx.*z - (1/fe)*At(A(z)) + 1/fe^2*(atr'*z)*atr;
nikcleju@2 240 # [dx, cgres, cgiter] = cgsolve(h11pfun, w1p, cgtol, cgmaxiter, 0);
nikcleju@2 241 # if (cgres > 1/2)
nikcleju@2 242 # disp('Cannot solve system. Returning previous iterate. (See Section 4 of notes for more information.)');
nikcleju@2 243 # xp = x; up = u;
nikcleju@2 244 # return
nikcleju@2 245 # end
nikcleju@2 246 # Adx = A(dx);
nikcleju@2 247 # else
nikcleju@2 248 # H11p = diag(sigx) - (1/fe)*AtA + (1/fe)^2*atr*atr';
nikcleju@2 249 # opts.POSDEF = true; opts.SYM = true;
nikcleju@2 250 # [dx,hcond] = linsolve(H11p, w1p, opts);
nikcleju@2 251 # if (hcond < 1e-14)
nikcleju@2 252 # disp('Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)');
nikcleju@2 253 # xp = x; up = u;
nikcleju@2 254 # return
nikcleju@2 255 # end
nikcleju@2 256 # Adx = A*dx;
nikcleju@2 257 # end
nikcleju@2 258 # du = (1./sig11).*ntgu - (sig12./sig11).*dx;
nikcleju@2 259 #
nikcleju@2 260 # # minimum step size that stays in the interior
nikcleju@2 261 # ifu1 = find((dx-du) > 0); ifu2 = find((-dx-du) > 0);
nikcleju@2 262 # aqe = Adx'*Adx; bqe = 2*r'*Adx; cqe = r'*r - epsilon^2;
nikcleju@2 263 # smax = min(1,min([...
nikcleju@2 264 # -fu1(ifu1)./(dx(ifu1)-du(ifu1)); -fu2(ifu2)./(-dx(ifu2)-du(ifu2)); ...
nikcleju@2 265 # (-bqe+sqrt(bqe^2-4*aqe*cqe))/(2*aqe)
nikcleju@2 266 # ]));
nikcleju@2 267 # s = (0.99)*smax;
nikcleju@2 268 #
nikcleju@2 269 # # backtracking line search
nikcleju@2 270 # suffdec = 0;
nikcleju@2 271 # backiter = 0;
nikcleju@2 272 # while (~suffdec)
nikcleju@2 273 # xp = x + s*dx; up = u + s*du; rp = r + s*Adx;
nikcleju@2 274 # fu1p = xp - up; fu2p = -xp - up; fep = 1/2*(rp'*rp - epsilon^2);
nikcleju@2 275 # fp = sum(up) - (1/tau)*(sum(log(-fu1p)) + sum(log(-fu2p)) + log(-fep));
nikcleju@2 276 # flin = f + alpha*s*(gradf'*[dx; du]);
nikcleju@2 277 # suffdec = (fp <= flin);
nikcleju@2 278 # s = beta*s;
nikcleju@2 279 # backiter = backiter + 1;
nikcleju@2 280 # if (backiter > 32)
nikcleju@2 281 # disp('Stuck on backtracking line search, returning previous iterate. (See Section 4 of notes for more information.)');
nikcleju@2 282 # xp = x; up = u;
nikcleju@2 283 # return
nikcleju@2 284 # end
nikcleju@2 285 # end
nikcleju@2 286 #
nikcleju@2 287 # # set up for next iteration
nikcleju@2 288 # x = xp; u = up; r = rp;
nikcleju@2 289 # fu1 = fu1p; fu2 = fu2p; fe = fep; f = fp;
nikcleju@2 290 #
nikcleju@2 291 # lambda2 = -(gradf'*[dx; du]);
nikcleju@2 292 # stepsize = s*norm([dx; du]);
nikcleju@2 293 # niter = niter + 1;
nikcleju@2 294 # done = (lambda2/2 < newtontol) | (niter >= newtonmaxiter);
nikcleju@2 295 #
nikcleju@2 296 # disp(sprintf('Newton iter = #d, Functional = #8.3f, Newton decrement = #8.3f, Stepsize = #8.3e', ...
nikcleju@2 297 # niter, f, lambda2/2, stepsize));
nikcleju@2 298 # if (largescale)
nikcleju@2 299 # disp(sprintf(' CG Res = #8.3e, CG Iter = #d', cgres, cgiter));
nikcleju@2 300 # else
nikcleju@2 301 # disp(sprintf(' H11p condition number = #8.3e', hcond));
nikcleju@2 302 # end
nikcleju@2 303 #
nikcleju@2 304 #end
nikcleju@2 305
nikcleju@2 306 # End of original Matab code
nikcleju@2 307 #----------------------------
nikcleju@2 308
nikcleju@2 309 # check if the matrix A is implicit or explicit
nikcleju@2 310 #largescale = isa(A,'function_handle');
nikcleju@2 311 if hasattr(A, '__call__'):
nikcleju@2 312 largescale = True
nikcleju@2 313 else:
nikcleju@2 314 largescale = False
nikcleju@2 315
nikcleju@2 316 # line search parameters
nikcleju@2 317 alpha = 0.01
nikcleju@2 318 beta = 0.5
nikcleju@2 319
nikcleju@2 320 #if (~largescale), AtA = A'*A; end
nikcleju@2 321 if not largescale:
nikcleju@2 322 AtA = np.dot(A.T,A)
nikcleju@2 323
nikcleju@2 324 # initial point
nikcleju@2 325 x = x0.copy()
nikcleju@2 326 u = u0.copy()
nikcleju@2 327 #if (largescale), r = A(x) - b; else r = A*x - b; end
nikcleju@2 328 if largescale:
nikcleju@2 329 r = A(x) - b
nikcleju@2 330 else:
nikcleju@2 331 r = np.dot(A,x) - b
nikcleju@2 332
nikcleju@2 333 fu1 = x - u
nikcleju@2 334 fu2 = -x - u
nikcleju@3 335 fe = 1.0/2*(np.vdot(r,r) - epsilon**2)
nikcleju@3 336 f = u.sum() - (1.0/tau)*(np.log(-fu1).sum() + np.log(-fu2).sum() + math.log(-fe))
nikcleju@2 337
nikcleju@2 338 niter = 0
nikcleju@2 339 done = 0
nikcleju@2 340 while not done:
nikcleju@2 341
nikcleju@2 342 #if (largescale), atr = At(r); else atr = A'*r; end
nikcleju@2 343 if largescale:
nikcleju@2 344 atr = At(r)
nikcleju@2 345 else:
nikcleju@2 346 atr = np.dot(A.T,r)
nikcleju@2 347
nikcleju@2 348 #ntgz = 1./fu1 - 1./fu2 + 1/fe*atr;
nikcleju@2 349 ntgz = 1.0/fu1 - 1.0/fu2 + 1.0/fe*atr
nikcleju@2 350 #ntgu = -tau - 1./fu1 - 1./fu2;
nikcleju@2 351 ntgu = -tau - 1.0/fu1 - 1.0/fu2
nikcleju@2 352 #gradf = -(1/tau)*[ntgz; ntgu];
nikcleju@2 353 gradf = -(1.0/tau)*np.concatenate((ntgz, ntgu),0)
nikcleju@2 354
nikcleju@2 355 #sig11 = 1./fu1.^2 + 1./fu2.^2;
nikcleju@2 356 sig11 = 1.0/(fu1**2) + 1.0/(fu2**2)
nikcleju@2 357 #sig12 = -1./fu1.^2 + 1./fu2.^2;
nikcleju@2 358 sig12 = -1.0/(fu1**2) + 1.0/(fu2**2)
nikcleju@2 359 #sigx = sig11 - sig12.^2./sig11;
nikcleju@2 360 sigx = sig11 - (sig12**2)/sig11
nikcleju@2 361
nikcleju@2 362 #w1p = ntgz - sig12./sig11.*ntgu;
nikcleju@2 363 w1p = ntgz - sig12/sig11*ntgu
nikcleju@2 364 if largescale:
nikcleju@2 365 #h11pfun = @(z) sigx.*z - (1/fe)*At(A(z)) + 1/fe^2*(atr'*z)*atr;
nikcleju@2 366 h11pfun = lambda z: sigx*z - (1.0/fe)*At(A(z)) + 1.0/(fe**2)*np.dot(np.dot(atr.T,z),atr)
nikcleju@2 367 dx,cgres,cgiter = cgsolve(h11pfun, w1p, cgtol, cgmaxiter, 0)
nikcleju@3 368 if (cgres > 1.0/2):
nikcleju@37 369 if verbose:
nikcleju@37 370 print 'Cannot solve system. Returning previous iterate. (See Section 4 of notes for more information.)'
nikcleju@2 371 xp = x.copy()
nikcleju@2 372 up = u.copy()
nikcleju@2 373 return xp,up,niter
nikcleju@2 374 #end
nikcleju@2 375 Adx = A(dx)
nikcleju@2 376 else:
nikcleju@2 377 #H11p = diag(sigx) - (1/fe)*AtA + (1/fe)^2*atr*atr';
nikcleju@3 378 # Attention: atr is column vector, so atr*atr' means outer(atr,atr)
nikcleju@3 379 H11p = np.diag(sigx) - (1.0/fe)*AtA + (1.0/fe)**2*np.outer(atr,atr)
nikcleju@2 380 #opts.POSDEF = true; opts.SYM = true;
nikcleju@2 381 #[dx,hcond] = linsolve(H11p, w1p, opts);
nikcleju@2 382 #if (hcond < 1e-14)
nikcleju@2 383 # disp('Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)');
nikcleju@2 384 # xp = x; up = u;
nikcleju@2 385 # return
nikcleju@2 386 #end
nikcleju@2 387 try:
nikcleju@2 388 dx = scipy.linalg.solve(H11p, w1p, sym_pos=True)
nikcleju@3 389 #dx = np.linalg.solve(H11p, w1p)
nikcleju@3 390 hcond = 1.0/np.linalg.cond(H11p)
nikcleju@2 391 except scipy.linalg.LinAlgError:
nikcleju@37 392 if verbose:
nikcleju@37 393 print 'Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)'
nikcleju@2 394 xp = x.copy()
nikcleju@2 395 up = u.copy()
nikcleju@2 396 return xp,up,niter
nikcleju@2 397 if hcond < 1e-14:
nikcleju@37 398 if verbose:
nikcleju@37 399 print 'Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)'
nikcleju@2 400 xp = x.copy()
nikcleju@2 401 up = u.copy()
nikcleju@2 402 return xp,up,niter
nikcleju@2 403
nikcleju@2 404 #Adx = A*dx;
nikcleju@2 405 Adx = np.dot(A,dx)
nikcleju@2 406 #end
nikcleju@2 407 #du = (1./sig11).*ntgu - (sig12./sig11).*dx;
nikcleju@2 408 du = (1.0/sig11)*ntgu - (sig12/sig11)*dx;
nikcleju@2 409
nikcleju@2 410 # minimum step size that stays in the interior
nikcleju@2 411 #ifu1 = find((dx-du) > 0); ifu2 = find((-dx-du) > 0);
nikcleju@2 412 ifu1 = np.nonzero((dx-du)>0)
nikcleju@2 413 ifu2 = np.nonzero((-dx-du)>0)
nikcleju@2 414 #aqe = Adx'*Adx; bqe = 2*r'*Adx; cqe = r'*r - epsilon^2;
nikcleju@2 415 aqe = np.dot(Adx.T,Adx)
nikcleju@2 416 bqe = 2*np.dot(r.T,Adx)
nikcleju@2 417 cqe = np.vdot(r,r) - epsilon**2
nikcleju@2 418 #smax = min(1,min([...
nikcleju@2 419 # -fu1(ifu1)./(dx(ifu1)-du(ifu1)); -fu2(ifu2)./(-dx(ifu2)-du(ifu2)); ...
nikcleju@2 420 # (-bqe+sqrt(bqe^2-4*aqe*cqe))/(2*aqe)
nikcleju@2 421 # ]));
nikcleju@3 422 smax = min(1,np.concatenate( (-fu1[ifu1]/(dx[ifu1]-du[ifu1]) , -fu2[ifu2]/(-dx[ifu2]-du[ifu2]) , np.array([ (-bqe + math.sqrt(bqe**2-4*aqe*cqe))/(2*aqe) ]) ) , 0).min())
nikcleju@3 423
nikcleju@3 424 s = 0.99 * smax
nikcleju@2 425
nikcleju@2 426 # backtracking line search
nikcleju@2 427 suffdec = 0
nikcleju@2 428 backiter = 0
nikcleju@2 429 while not suffdec:
nikcleju@2 430 #xp = x + s*dx; up = u + s*du; rp = r + s*Adx;
nikcleju@2 431 xp = x + s*dx
nikcleju@2 432 up = u + s*du
nikcleju@2 433 rp = r + s*Adx
nikcleju@2 434 #fu1p = xp - up; fu2p = -xp - up; fep = 1/2*(rp'*rp - epsilon^2);
nikcleju@2 435 fu1p = xp - up
nikcleju@2 436 fu2p = -xp - up
nikcleju@3 437 fep = 0.5*(np.vdot(rp,rp) - epsilon**2)
nikcleju@2 438 #fp = sum(up) - (1/tau)*(sum(log(-fu1p)) + sum(log(-fu2p)) + log(-fep));
nikcleju@3 439 fp = up.sum() - (1.0/tau)*(np.log(-fu1p).sum() + np.log(-fu2p).sum() + math.log(-fep))
nikcleju@2 440 #flin = f + alpha*s*(gradf'*[dx; du]);
nikcleju@2 441 flin = f + alpha*s*np.dot(gradf.T , np.concatenate((dx,du),0))
nikcleju@2 442 #suffdec = (fp <= flin);
nikcleju@2 443 if fp <= flin:
nikcleju@2 444 suffdec = True
nikcleju@2 445 else:
nikcleju@2 446 suffdec = False
nikcleju@2 447
nikcleju@2 448 s = beta*s
nikcleju@2 449 backiter = backiter + 1
nikcleju@2 450 if (backiter > 32):
nikcleju@37 451 if verbose:
nikcleju@37 452 print 'Stuck on backtracking line search, returning previous iterate. (See Section 4 of notes for more information.)'
nikcleju@2 453 xp = x.copy()
nikcleju@2 454 up = u.copy()
nikcleju@2 455 return xp,up,niter
nikcleju@2 456 #end
nikcleju@2 457 #end
nikcleju@2 458
nikcleju@2 459 # set up for next iteration
nikcleju@2 460 #x = xp; u = up; r = rp;
nikcleju@2 461 x = xp.copy()
nikcleju@2 462 u = up.copy()
nikcleju@2 463 r = rp.copy()
nikcleju@2 464 #fu1 = fu1p; fu2 = fu2p; fe = fep; f = fp;
nikcleju@2 465 fu1 = fu1p.copy()
nikcleju@2 466 fu2 = fu2p.copy()
nikcleju@2 467 fe = fep
nikcleju@2 468 f = fp
nikcleju@2 469
nikcleju@2 470 #lambda2 = -(gradf'*[dx; du]);
nikcleju@3 471 lambda2 = -np.dot(gradf.T , np.concatenate((dx,du),0))
nikcleju@2 472 #stepsize = s*norm([dx; du]);
nikcleju@2 473 stepsize = s * np.linalg.norm(np.concatenate((dx,du),0))
nikcleju@2 474 niter = niter + 1
nikcleju@2 475 #done = (lambda2/2 < newtontol) | (niter >= newtonmaxiter);
nikcleju@2 476 if lambda2/2.0 < newtontol or niter >= newtonmaxiter:
nikcleju@2 477 done = 1
nikcleju@2 478 else:
nikcleju@2 479 done = 0
nikcleju@2 480
nikcleju@2 481 #disp(sprintf('Newton iter = #d, Functional = #8.3f, Newton decrement = #8.3f, Stepsize = #8.3e', ...
nikcleju@37 482 if verbose:
nikcleju@37 483 print 'Newton iter = ',niter,', Functional = ',f,', Newton decrement = ',lambda2/2.0,', Stepsize = ',stepsize
nikcleju@2 484
nikcleju@37 485 if verbose:
nikcleju@37 486 if largescale:
nikcleju@37 487 print ' CG Res = ',cgres,', CG Iter = ',cgiter
nikcleju@37 488 else:
nikcleju@37 489 print ' H11p condition number = ',hcond
nikcleju@2 490 #end
nikcleju@2 491
nikcleju@2 492 #end
nikcleju@2 493 return xp,up,niter
nikcleju@2 494
nikcleju@2 495 #function xp = l1qc_logbarrier(x0, A, At, b, epsilon, lbtol, mu, cgtol, cgmaxiter)
nikcleju@37 496 def l1qc_logbarrier(x0, A, At, b, epsilon, lbtol=1e-3, mu=10, cgtol=1e-8, cgmaxiter=200, verbose=False):
nikcleju@2 497 # Solve quadratically constrained l1 minimization:
nikcleju@2 498 # min ||x||_1 s.t. ||Ax - b||_2 <= \epsilon
nikcleju@2 499 #
nikcleju@2 500 # Reformulate as the second-order cone program
nikcleju@2 501 # min_{x,u} sum(u) s.t. x - u <= 0,
nikcleju@2 502 # -x - u <= 0,
nikcleju@2 503 # 1/2(||Ax-b||^2 - \epsilon^2) <= 0
nikcleju@2 504 # and use a log barrier algorithm.
nikcleju@2 505 #
nikcleju@2 506 # Usage: xp = l1qc_logbarrier(x0, A, At, b, epsilon, lbtol, mu, cgtol, cgmaxiter)
nikcleju@2 507 #
nikcleju@2 508 # x0 - Nx1 vector, initial point.
nikcleju@2 509 #
nikcleju@2 510 # A - Either a handle to a function that takes a N vector and returns a K
nikcleju@2 511 # vector , or a KxN matrix. If A is a function handle, the algorithm
nikcleju@2 512 # operates in "largescale" mode, solving the Newton systems via the
nikcleju@2 513 # Conjugate Gradients algorithm.
nikcleju@2 514 #
nikcleju@2 515 # At - Handle to a function that takes a K vector and returns an N vector.
nikcleju@2 516 # If A is a KxN matrix, At is ignored.
nikcleju@2 517 #
nikcleju@2 518 # b - Kx1 vector of observations.
nikcleju@2 519 #
nikcleju@2 520 # epsilon - scalar, constraint relaxation parameter
nikcleju@2 521 #
nikcleju@2 522 # lbtol - The log barrier algorithm terminates when the duality gap <= lbtol.
nikcleju@2 523 # Also, the number of log barrier iterations is completely
nikcleju@2 524 # determined by lbtol.
nikcleju@2 525 # Default = 1e-3.
nikcleju@2 526 #
nikcleju@2 527 # mu - Factor by which to increase the barrier constant at each iteration.
nikcleju@2 528 # Default = 10.
nikcleju@2 529 #
nikcleju@2 530 # cgtol - Tolerance for Conjugate Gradients; ignored if A is a matrix.
nikcleju@2 531 # Default = 1e-8.
nikcleju@2 532 #
nikcleju@2 533 # cgmaxiter - Maximum number of iterations for Conjugate Gradients; ignored
nikcleju@2 534 # if A is a matrix.
nikcleju@2 535 # Default = 200.
nikcleju@2 536 #
nikcleju@2 537 # Written by: Justin Romberg, Caltech
nikcleju@2 538 # Email: jrom@acm.caltech.edu
nikcleju@2 539 # Created: October 2005
nikcleju@2 540 #
nikcleju@2 541
nikcleju@2 542 #---------------------
nikcleju@2 543 # Original Matab code:
nikcleju@2 544
nikcleju@2 545 #largescale = isa(A,'function_handle');
nikcleju@2 546 #
nikcleju@2 547 #if (nargin < 6), lbtol = 1e-3; end
nikcleju@2 548 #if (nargin < 7), mu = 10; end
nikcleju@2 549 #if (nargin < 8), cgtol = 1e-8; end
nikcleju@2 550 #if (nargin < 9), cgmaxiter = 200; end
nikcleju@2 551 #
nikcleju@2 552 #newtontol = lbtol;
nikcleju@2 553 #newtonmaxiter = 50;
nikcleju@2 554 #
nikcleju@2 555 #N = length(x0);
nikcleju@2 556 #
nikcleju@2 557 ## starting point --- make sure that it is feasible
nikcleju@2 558 #if (largescale)
nikcleju@2 559 # if (norm(A(x0)-b) > epsilon)
nikcleju@2 560 # disp('Starting point infeasible; using x0 = At*inv(AAt)*y.');
nikcleju@2 561 # AAt = @(z) A(At(z));
nikcleju@2 562 # w = cgsolve(AAt, b, cgtol, cgmaxiter, 0);
nikcleju@2 563 # if (cgres > 1/2)
nikcleju@2 564 # disp('A*At is ill-conditioned: cannot find starting point');
nikcleju@2 565 # xp = x0;
nikcleju@2 566 # return;
nikcleju@2 567 # end
nikcleju@2 568 # x0 = At(w);
nikcleju@2 569 # end
nikcleju@2 570 #else
nikcleju@2 571 # if (norm(A*x0-b) > epsilon)
nikcleju@2 572 # disp('Starting point infeasible; using x0 = At*inv(AAt)*y.');
nikcleju@2 573 # opts.POSDEF = true; opts.SYM = true;
nikcleju@2 574 # [w, hcond] = linsolve(A*A', b, opts);
nikcleju@2 575 # if (hcond < 1e-14)
nikcleju@2 576 # disp('A*At is ill-conditioned: cannot find starting point');
nikcleju@2 577 # xp = x0;
nikcleju@2 578 # return;
nikcleju@2 579 # end
nikcleju@2 580 # x0 = A'*w;
nikcleju@2 581 # end
nikcleju@2 582 #end
nikcleju@2 583 #x = x0;
nikcleju@2 584 #u = (0.95)*abs(x0) + (0.10)*max(abs(x0));
nikcleju@2 585 #
nikcleju@2 586 #disp(sprintf('Original l1 norm = #.3f, original functional = #.3f', sum(abs(x0)), sum(u)));
nikcleju@2 587 #
nikcleju@2 588 ## choose initial value of tau so that the duality gap after the first
nikcleju@2 589 ## step will be about the origial norm
nikcleju@2 590 #tau = max((2*N+1)/sum(abs(x0)), 1);
nikcleju@2 591 #
nikcleju@2 592 #lbiter = ceil((log(2*N+1)-log(lbtol)-log(tau))/log(mu));
nikcleju@2 593 #disp(sprintf('Number of log barrier iterations = #d\n', lbiter));
nikcleju@2 594 #
nikcleju@2 595 #totaliter = 0;
nikcleju@2 596 #
nikcleju@2 597 ## Added by Nic
nikcleju@2 598 #if lbiter == 0
nikcleju@2 599 # xp = zeros(size(x0));
nikcleju@2 600 #end
nikcleju@2 601 #
nikcleju@2 602 #for ii = 1:lbiter
nikcleju@2 603 #
nikcleju@2 604 # [xp, up, ntiter] = l1qc_newton(x, u, A, At, b, epsilon, tau, newtontol, newtonmaxiter, cgtol, cgmaxiter);
nikcleju@2 605 # totaliter = totaliter + ntiter;
nikcleju@2 606 #
nikcleju@2 607 # disp(sprintf('\nLog barrier iter = #d, l1 = #.3f, functional = #8.3f, tau = #8.3e, total newton iter = #d\n', ...
nikcleju@2 608 # ii, sum(abs(xp)), sum(up), tau, totaliter));
nikcleju@2 609 #
nikcleju@2 610 # x = xp;
nikcleju@2 611 # u = up;
nikcleju@2 612 #
nikcleju@2 613 # tau = mu*tau;
nikcleju@2 614 #
nikcleju@2 615 #end
nikcleju@2 616 #
nikcleju@2 617 # End of original Matab code
nikcleju@2 618 #----------------------------
nikcleju@2 619
nikcleju@37 620 # Check if epsilon > 0. If epsilon is 0, the algorithm fails. You should run the algo with equality constraint instead
nikcleju@37 621 if epsilon == 0:
nikcleju@37 622 raise l1qcInputValueError('Epsilon should be > 0!')
nikcleju@37 623
nikcleju@2 624 #largescale = isa(A,'function_handle');
nikcleju@2 625 if hasattr(A, '__call__'):
nikcleju@2 626 largescale = True
nikcleju@2 627 else:
nikcleju@2 628 largescale = False
nikcleju@2 629
nikcleju@2 630 # if (nargin < 6), lbtol = 1e-3; end
nikcleju@2 631 # if (nargin < 7), mu = 10; end
nikcleju@2 632 # if (nargin < 8), cgtol = 1e-8; end
nikcleju@2 633 # if (nargin < 9), cgmaxiter = 200; end
nikcleju@2 634 # Nic: added them as optional parameteres
nikcleju@2 635
nikcleju@2 636 newtontol = lbtol
nikcleju@2 637 newtonmaxiter = 50
nikcleju@2 638
nikcleju@2 639 #N = length(x0);
nikcleju@3 640 N = x0.size
nikcleju@2 641
nikcleju@2 642 # starting point --- make sure that it is feasible
nikcleju@2 643 if largescale:
nikcleju@2 644 if np.linalg.norm(A(x0) - b) > epsilon:
nikcleju@37 645 if verbose:
nikcleju@37 646 print 'Starting point infeasible; using x0 = At*inv(AAt)*y.'
nikcleju@2 647 #AAt = @(z) A(At(z));
nikcleju@2 648 AAt = lambda z: A(At(z))
nikcleju@2 649 # TODO: implement cgsolve
nikcleju@2 650 w,cgres,cgiter = cgsolve(AAt, b, cgtol, cgmaxiter, 0)
nikcleju@3 651 if (cgres > 1.0/2):
nikcleju@37 652 if verbose:
nikcleju@37 653 print 'A*At is ill-conditioned: cannot find starting point'
nikcleju@2 654 xp = x0.copy()
nikcleju@2 655 return xp
nikcleju@2 656 #end
nikcleju@2 657 x0 = At(w)
nikcleju@2 658 #end
nikcleju@2 659 else:
nikcleju@2 660 if np.linalg.norm( np.dot(A,x0) - b ) > epsilon:
nikcleju@37 661 if verbose:
nikcleju@37 662 print 'Starting point infeasible; using x0 = At*inv(AAt)*y.'
nikcleju@2 663 #opts.POSDEF = true; opts.SYM = true;
nikcleju@2 664 #[w, hcond] = linsolve(A*A', b, opts);
nikcleju@2 665 #if (hcond < 1e-14)
nikcleju@2 666 # disp('A*At is ill-conditioned: cannot find starting point');
nikcleju@2 667 # xp = x0;
nikcleju@2 668 # return;
nikcleju@2 669 #end
nikcleju@2 670 try:
nikcleju@2 671 w = scipy.linalg.solve(np.dot(A,A.T), b, sym_pos=True)
nikcleju@3 672 #w = np.linalg.solve(np.dot(A,A.T), b)
nikcleju@37 673 hcond = 1.0/np.linalg.cond(np.dot(A,A.T))
nikcleju@2 674 except scipy.linalg.LinAlgError:
nikcleju@37 675 if verbose:
nikcleju@37 676 print 'A*At is ill-conditioned: cannot find starting point'
nikcleju@2 677 xp = x0.copy()
nikcleju@2 678 return xp
nikcleju@2 679 if hcond < 1e-14:
nikcleju@37 680 if verbose:
nikcleju@37 681 print 'A*At is ill-conditioned: cannot find starting point'
nikcleju@2 682 xp = x0.copy()
nikcleju@2 683 return xp
nikcleju@2 684 #x0 = A'*w;
nikcleju@2 685 x0 = np.dot(A.T, w)
nikcleju@2 686 #end
nikcleju@2 687 #end
nikcleju@2 688 x = x0.copy()
nikcleju@2 689 u = (0.95)*np.abs(x0) + (0.10)*np.abs(x0).max()
nikcleju@2 690
nikcleju@2 691 #disp(sprintf('Original l1 norm = #.3f, original functional = #.3f', sum(abs(x0)), sum(u)));
nikcleju@37 692 if verbose:
nikcleju@37 693 print 'Original l1 norm = ',np.abs(x0).sum(),'original functional = ',u.sum()
nikcleju@2 694
nikcleju@2 695 # choose initial value of tau so that the duality gap after the first
nikcleju@2 696 # step will be about the origial norm
nikcleju@3 697 tau = max(((2*N+1.0)/np.abs(x0).sum()), 1)
nikcleju@2 698
nikcleju@2 699 lbiter = math.ceil((math.log(2*N+1)-math.log(lbtol)-math.log(tau))/math.log(mu))
nikcleju@2 700 #disp(sprintf('Number of log barrier iterations = #d\n', lbiter));
nikcleju@37 701 if verbose:
nikcleju@37 702 print 'Number of log barrier iterations = ',lbiter
nikcleju@2 703
nikcleju@2 704 totaliter = 0
nikcleju@2 705
nikcleju@2 706 # Added by Nic, to fix some crashing
nikcleju@2 707 if lbiter == 0:
nikcleju@2 708 xp = np.zeros(x0.size)
nikcleju@2 709 #end
nikcleju@2 710
nikcleju@2 711 #for ii = 1:lbiter
nikcleju@2 712 for ii in np.arange(lbiter):
nikcleju@2 713
nikcleju@2 714 xp,up,ntiter = l1qc_newton(x, u, A, At, b, epsilon, tau, newtontol, newtonmaxiter, cgtol, cgmaxiter)
nikcleju@2 715 totaliter = totaliter + ntiter
nikcleju@2 716
nikcleju@2 717 #disp(sprintf('\nLog barrier iter = #d, l1 = #.3f, functional = #8.3f, tau = #8.3e, total newton iter = #d\n', ...
nikcleju@2 718 # ii, sum(abs(xp)), sum(up), tau, totaliter));
nikcleju@37 719 if verbose:
nikcleju@37 720 print 'Log barrier iter = ',ii,', l1 = ',np.abs(xp).sum(),', functional = ',up.sum(),', tau = ',tau,', total newton iter = ',totaliter
nikcleju@2 721 x = xp.copy()
nikcleju@2 722 u = up.copy()
nikcleju@2 723
nikcleju@2 724 tau = mu*tau
nikcleju@2 725
nikcleju@2 726 #end
nikcleju@2 727 return xp
nikcleju@2 728