nikcleju@2
|
1 % l1qc_newton.m
|
nikcleju@2
|
2 %
|
nikcleju@2
|
3 % Newton algorithm for log-barrier subproblems for l1 minimization
|
nikcleju@2
|
4 % with quadratic constraints.
|
nikcleju@2
|
5 %
|
nikcleju@2
|
6 % Usage:
|
nikcleju@2
|
7 % [xp,up,niter] = l1qc_newton(x0, u0, A, At, b, epsilon, tau,
|
nikcleju@2
|
8 % newtontol, newtonmaxiter, cgtol, cgmaxiter)
|
nikcleju@2
|
9 %
|
nikcleju@2
|
10 % x0,u0 - starting points
|
nikcleju@2
|
11 %
|
nikcleju@2
|
12 % A - Either a handle to a function that takes a N vector and returns a K
|
nikcleju@2
|
13 % vector , or a KxN matrix. If A is a function handle, the algorithm
|
nikcleju@2
|
14 % operates in "largescale" mode, solving the Newton systems via the
|
nikcleju@2
|
15 % Conjugate Gradients algorithm.
|
nikcleju@2
|
16 %
|
nikcleju@2
|
17 % At - Handle to a function that takes a K vector and returns an N vector.
|
nikcleju@2
|
18 % If A is a KxN matrix, At is ignored.
|
nikcleju@2
|
19 %
|
nikcleju@2
|
20 % b - Kx1 vector of observations.
|
nikcleju@2
|
21 %
|
nikcleju@2
|
22 % epsilon - scalar, constraint relaxation parameter
|
nikcleju@2
|
23 %
|
nikcleju@2
|
24 % tau - Log barrier parameter.
|
nikcleju@2
|
25 %
|
nikcleju@2
|
26 % newtontol - Terminate when the Newton decrement is <= newtontol.
|
nikcleju@2
|
27 % Default = 1e-3.
|
nikcleju@2
|
28 %
|
nikcleju@2
|
29 % newtonmaxiter - Maximum number of iterations.
|
nikcleju@2
|
30 % Default = 50.
|
nikcleju@2
|
31 %
|
nikcleju@2
|
32 % cgtol - Tolerance for Conjugate Gradients; ignored if A is a matrix.
|
nikcleju@2
|
33 % Default = 1e-8.
|
nikcleju@2
|
34 %
|
nikcleju@2
|
35 % cgmaxiter - Maximum number of iterations for Conjugate Gradients; ignored
|
nikcleju@2
|
36 % if A is a matrix.
|
nikcleju@2
|
37 % Default = 200.
|
nikcleju@2
|
38 %
|
nikcleju@2
|
39 % Written by: Justin Romberg, Caltech
|
nikcleju@2
|
40 % Email: jrom@acm.caltech.edu
|
nikcleju@2
|
41 % Created: October 2005
|
nikcleju@2
|
42 %
|
nikcleju@2
|
43
|
nikcleju@2
|
44
|
nikcleju@2
|
45 function [xp, up, niter] = l1qc_newton(x0, u0, A, At, b, epsilon, tau, newtontol, newtonmaxiter, cgtol, cgmaxiter)
|
nikcleju@2
|
46
|
nikcleju@2
|
47 % check if the matrix A is implicit or explicit
|
nikcleju@2
|
48 largescale = isa(A,'function_handle');
|
nikcleju@2
|
49
|
nikcleju@2
|
50 % line search parameters
|
nikcleju@2
|
51 alpha = 0.01;
|
nikcleju@2
|
52 beta = 0.5;
|
nikcleju@2
|
53
|
nikcleju@2
|
54 if (~largescale), AtA = A'*A; end
|
nikcleju@2
|
55
|
nikcleju@2
|
56 % initial point
|
nikcleju@2
|
57 x = x0;
|
nikcleju@2
|
58 u = u0;
|
nikcleju@2
|
59 if (largescale), r = A(x) - b; else r = A*x - b; end
|
nikcleju@2
|
60 fu1 = x - u;
|
nikcleju@2
|
61 fu2 = -x - u;
|
nikcleju@2
|
62 fe = 1/2*(r'*r - epsilon^2);
|
nikcleju@2
|
63 f = sum(u) - (1/tau)*(sum(log(-fu1)) + sum(log(-fu2)) + log(-fe));
|
nikcleju@2
|
64
|
nikcleju@2
|
65 niter = 0;
|
nikcleju@2
|
66 done = 0;
|
nikcleju@2
|
67 while (~done)
|
nikcleju@2
|
68
|
nikcleju@2
|
69 if (largescale), atr = At(r); else atr = A'*r; end
|
nikcleju@2
|
70
|
nikcleju@2
|
71 ntgz = 1./fu1 - 1./fu2 + 1/fe*atr;
|
nikcleju@2
|
72 ntgu = -tau - 1./fu1 - 1./fu2;
|
nikcleju@2
|
73 gradf = -(1/tau)*[ntgz; ntgu];
|
nikcleju@2
|
74
|
nikcleju@2
|
75 sig11 = 1./fu1.^2 + 1./fu2.^2;
|
nikcleju@2
|
76 sig12 = -1./fu1.^2 + 1./fu2.^2;
|
nikcleju@2
|
77 sigx = sig11 - sig12.^2./sig11;
|
nikcleju@2
|
78
|
nikcleju@2
|
79 w1p = ntgz - sig12./sig11.*ntgu;
|
nikcleju@2
|
80 if (largescale)
|
nikcleju@2
|
81 h11pfun = @(z) sigx.*z - (1/fe)*At(A(z)) + 1/fe^2*(atr'*z)*atr;
|
nikcleju@2
|
82 [dx, cgres, cgiter] = cgsolve(h11pfun, w1p, cgtol, cgmaxiter, 0);
|
nikcleju@2
|
83 if (cgres > 1/2)
|
nikcleju@2
|
84 disp('Cannot solve system. Returning previous iterate. (See Section 4 of notes for more information.)');
|
nikcleju@2
|
85 xp = x; up = u;
|
nikcleju@2
|
86 return
|
nikcleju@2
|
87 end
|
nikcleju@2
|
88 Adx = A(dx);
|
nikcleju@2
|
89 else
|
nikcleju@2
|
90 H11p = diag(sigx) - (1/fe)*AtA + (1/fe)^2*atr*atr';
|
nikcleju@2
|
91 opts.POSDEF = true; opts.SYM = true;
|
nikcleju@2
|
92 [dx,hcond] = linsolve(H11p, w1p, opts);
|
nikcleju@2
|
93 if (hcond < 1e-14)
|
nikcleju@2
|
94 disp('Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)');
|
nikcleju@2
|
95 xp = x; up = u;
|
nikcleju@2
|
96 return
|
nikcleju@2
|
97 end
|
nikcleju@2
|
98 Adx = A*dx;
|
nikcleju@2
|
99 end
|
nikcleju@2
|
100 du = (1./sig11).*ntgu - (sig12./sig11).*dx;
|
nikcleju@2
|
101
|
nikcleju@2
|
102 % minimum step size that stays in the interior
|
nikcleju@2
|
103 ifu1 = find((dx-du) > 0); ifu2 = find((-dx-du) > 0);
|
nikcleju@2
|
104 aqe = Adx'*Adx; bqe = 2*r'*Adx; cqe = r'*r - epsilon^2;
|
nikcleju@2
|
105 smax = min(1,min([...
|
nikcleju@2
|
106 -fu1(ifu1)./(dx(ifu1)-du(ifu1)); -fu2(ifu2)./(-dx(ifu2)-du(ifu2)); ...
|
nikcleju@2
|
107 (-bqe+sqrt(bqe^2-4*aqe*cqe))/(2*aqe)
|
nikcleju@2
|
108 ]));
|
nikcleju@2
|
109 s = (0.99)*smax;
|
nikcleju@2
|
110
|
nikcleju@2
|
111 % backtracking line search
|
nikcleju@2
|
112 suffdec = 0;
|
nikcleju@2
|
113 backiter = 0;
|
nikcleju@2
|
114 while (~suffdec)
|
nikcleju@2
|
115 xp = x + s*dx; up = u + s*du; rp = r + s*Adx;
|
nikcleju@2
|
116 fu1p = xp - up; fu2p = -xp - up; fep = 1/2*(rp'*rp - epsilon^2);
|
nikcleju@2
|
117 fp = sum(up) - (1/tau)*(sum(log(-fu1p)) + sum(log(-fu2p)) + log(-fep));
|
nikcleju@2
|
118 flin = f + alpha*s*(gradf'*[dx; du]);
|
nikcleju@2
|
119 suffdec = (fp <= flin);
|
nikcleju@2
|
120 s = beta*s;
|
nikcleju@2
|
121 backiter = backiter + 1;
|
nikcleju@2
|
122 if (backiter > 32)
|
nikcleju@2
|
123 disp('Stuck on backtracking line search, returning previous iterate. (See Section 4 of notes for more information.)');
|
nikcleju@2
|
124 xp = x; up = u;
|
nikcleju@2
|
125 return
|
nikcleju@2
|
126 end
|
nikcleju@2
|
127 end
|
nikcleju@2
|
128
|
nikcleju@2
|
129 % set up for next iteration
|
nikcleju@2
|
130 x = xp; u = up; r = rp;
|
nikcleju@2
|
131 fu1 = fu1p; fu2 = fu2p; fe = fep; f = fp;
|
nikcleju@2
|
132
|
nikcleju@2
|
133 lambda2 = -(gradf'*[dx; du]);
|
nikcleju@2
|
134 stepsize = s*norm([dx; du]);
|
nikcleju@2
|
135 niter = niter + 1;
|
nikcleju@2
|
136 done = (lambda2/2 < newtontol) | (niter >= newtonmaxiter);
|
nikcleju@2
|
137
|
nikcleju@2
|
138 disp(sprintf('Newton iter = %d, Functional = %8.3f, Newton decrement = %8.3f, Stepsize = %8.3e', ...
|
nikcleju@2
|
139 niter, f, lambda2/2, stepsize));
|
nikcleju@2
|
140 if (largescale)
|
nikcleju@2
|
141 disp(sprintf(' CG Res = %8.3e, CG Iter = %d', cgres, cgiter));
|
nikcleju@2
|
142 else
|
nikcleju@2
|
143 disp(sprintf(' H11p condition number = %8.3e', hcond));
|
nikcleju@2
|
144 end
|
nikcleju@2
|
145
|
nikcleju@2
|
146 end
|
nikcleju@2
|
147
|
nikcleju@2
|
148
|
nikcleju@2
|
149
|
nikcleju@2
|
150
|