annotate matlab/BP/l1eq_pd.m @ 68:cab8a215f9a1 tip

Minor
author Nic Cleju <nikcleju@gmail.com>
date Tue, 09 Jul 2013 14:50:09 +0300
parents e684f76c1969
children
rev   line source
nikcleju@62 1 % l1eq_pd.m
nikcleju@62 2 %
nikcleju@62 3 % Solve
nikcleju@62 4 % min_x ||x||_1 s.t. Ax = b
nikcleju@62 5 %
nikcleju@62 6 % Recast as linear program
nikcleju@62 7 % min_{x,u} sum(u) s.t. -u <= x <= u, Ax=b
nikcleju@62 8 % and use primal-dual interior point method
nikcleju@62 9 %
nikcleju@62 10 % Usage: xp = l1eq_pd(x0, A, At, b, pdtol, pdmaxiter, cgtol, cgmaxiter)
nikcleju@62 11 %
nikcleju@62 12 % x0 - Nx1 vector, initial point.
nikcleju@62 13 %
nikcleju@62 14 % A - Either a handle to a function that takes a N vector and returns a K
nikcleju@62 15 % vector , or a KxN matrix. If A is a function handle, the algorithm
nikcleju@62 16 % operates in "largescale" mode, solving the Newton systems via the
nikcleju@62 17 % Conjugate Gradients algorithm.
nikcleju@62 18 %
nikcleju@62 19 % At - Handle to a function that takes a K vector and returns an N vector.
nikcleju@62 20 % If A is a KxN matrix, At is ignored.
nikcleju@62 21 %
nikcleju@62 22 % b - Kx1 vector of observations.
nikcleju@62 23 %
nikcleju@62 24 % pdtol - Tolerance for primal-dual algorithm (algorithm terminates if
nikcleju@62 25 % the duality gap is less than pdtol).
nikcleju@62 26 % Default = 1e-3.
nikcleju@62 27 %
nikcleju@62 28 % pdmaxiter - Maximum number of primal-dual iterations.
nikcleju@62 29 % Default = 50.
nikcleju@62 30 %
nikcleju@62 31 % cgtol - Tolerance for Conjugate Gradients; ignored if A is a matrix.
nikcleju@62 32 % Default = 1e-8.
nikcleju@62 33 %
nikcleju@62 34 % cgmaxiter - Maximum number of iterations for Conjugate Gradients; ignored
nikcleju@62 35 % if A is a matrix.
nikcleju@62 36 % Default = 200.
nikcleju@62 37 %
nikcleju@62 38 % Written by: Justin Romberg, Caltech
nikcleju@62 39 % Email: jrom@acm.caltech.edu
nikcleju@62 40 % Created: October 2005
nikcleju@62 41 %
nikcleju@62 42
nikcleju@62 43 function xp = l1eq_pd(x0, A, At, b, pdtol, pdmaxiter, cgtol, cgmaxiter)
nikcleju@62 44
nikcleju@62 45 largescale = isa(A,'function_handle');
nikcleju@62 46
nikcleju@62 47 if (nargin < 5), pdtol = 1e-3; end
nikcleju@62 48 if (nargin < 6), pdmaxiter = 50; end
nikcleju@62 49 if (nargin < 7), cgtol = 1e-8; end
nikcleju@62 50 if (nargin < 8), cgmaxiter = 200; end
nikcleju@62 51
nikcleju@62 52 N = length(x0);
nikcleju@62 53
nikcleju@62 54 alpha = 0.01;
nikcleju@62 55 beta = 0.5;
nikcleju@62 56 mu = 10;
nikcleju@62 57
nikcleju@62 58 gradf0 = [zeros(N,1); ones(N,1)];
nikcleju@62 59
nikcleju@62 60 % starting point --- make sure that it is feasible
nikcleju@62 61 if (largescale)
nikcleju@62 62 if (norm(A(x0)-b)/norm(b) > cgtol)
nikcleju@62 63 disp('Starting point infeasible; using x0 = At*inv(AAt)*y.');
nikcleju@62 64 AAt = @(z) A(At(z));
nikcleju@62 65 [w, cgres, cgiter] = cgsolve(AAt, b, cgtol, cgmaxiter, 0);
nikcleju@62 66 if (cgres > 1/2)
nikcleju@62 67 disp('A*At is ill-conditioned: cannot find starting point');
nikcleju@62 68 xp = x0;
nikcleju@62 69 return;
nikcleju@62 70 end
nikcleju@62 71 x0 = At(w);
nikcleju@62 72 end
nikcleju@62 73 else
nikcleju@62 74 if (norm(A*x0-b)/norm(b) > cgtol)
nikcleju@62 75 disp('Starting point infeasible; using x0 = At*inv(AAt)*y.');
nikcleju@62 76 opts.POSDEF = true; opts.SYM = true;
nikcleju@62 77 [w, hcond] = linsolve(A*A', b, opts);
nikcleju@62 78 if (hcond < 1e-14)
nikcleju@62 79 disp('A*At is ill-conditioned: cannot find starting point');
nikcleju@62 80 xp = x0;
nikcleju@62 81 return;
nikcleju@62 82 end
nikcleju@62 83 x0 = A'*w;
nikcleju@62 84 end
nikcleju@62 85 end
nikcleju@62 86 x = x0;
nikcleju@62 87 u = (0.95)*abs(x0) + (0.10)*max(abs(x0));
nikcleju@62 88
nikcleju@62 89 % set up for the first iteration
nikcleju@62 90 fu1 = x - u;
nikcleju@62 91 fu2 = -x - u;
nikcleju@62 92 lamu1 = -1./fu1;
nikcleju@62 93 lamu2 = -1./fu2;
nikcleju@62 94 if (largescale)
nikcleju@62 95 v = -A(lamu1-lamu2);
nikcleju@62 96 Atv = At(v);
nikcleju@62 97 rpri = A(x) - b;
nikcleju@62 98 else
nikcleju@62 99 v = -A*(lamu1-lamu2);
nikcleju@62 100 Atv = A'*v;
nikcleju@62 101 rpri = A*x - b;
nikcleju@62 102 end
nikcleju@62 103
nikcleju@62 104 sdg = -(fu1'*lamu1 + fu2'*lamu2);
nikcleju@62 105 tau = mu*2*N/sdg;
nikcleju@62 106
nikcleju@62 107 rcent = [-lamu1.*fu1; -lamu2.*fu2] - (1/tau);
nikcleju@62 108 rdual = gradf0 + [lamu1-lamu2; -lamu1-lamu2] + [Atv; zeros(N,1)];
nikcleju@62 109 resnorm = norm([rdual; rcent; rpri]);
nikcleju@62 110
nikcleju@62 111 pditer = 0;
nikcleju@62 112 done = (sdg < pdtol) | (pditer >= pdmaxiter);
nikcleju@62 113 while (~done)
nikcleju@62 114
nikcleju@62 115 pditer = pditer + 1;
nikcleju@62 116
nikcleju@62 117 w1 = -1/tau*(-1./fu1 + 1./fu2) - Atv;
nikcleju@62 118 w2 = -1 - 1/tau*(1./fu1 + 1./fu2);
nikcleju@62 119 w3 = -rpri;
nikcleju@62 120
nikcleju@62 121 sig1 = -lamu1./fu1 - lamu2./fu2;
nikcleju@62 122 sig2 = lamu1./fu1 - lamu2./fu2;
nikcleju@62 123 sigx = sig1 - sig2.^2./sig1;
nikcleju@62 124
nikcleju@62 125 if (largescale)
nikcleju@62 126 w1p = w3 - A(w1./sigx - w2.*sig2./(sigx.*sig1));
nikcleju@62 127 h11pfun = @(z) -A(1./sigx.*At(z));
nikcleju@62 128 [dv, cgres, cgiter] = cgsolve(h11pfun, w1p, cgtol, cgmaxiter, 0);
nikcleju@62 129 if (cgres > 1/2)
nikcleju@62 130 disp('Cannot solve system. Returning previous iterate. (See Section 4 of notes for more information.)');
nikcleju@62 131 xp = x;
nikcleju@62 132 return
nikcleju@62 133 end
nikcleju@62 134 dx = (w1 - w2.*sig2./sig1 - At(dv))./sigx;
nikcleju@62 135 Adx = A(dx);
nikcleju@62 136 Atdv = At(dv);
nikcleju@62 137 else
nikcleju@62 138 w1p = -(w3 - A*(w1./sigx - w2.*sig2./(sigx.*sig1)));
nikcleju@62 139 H11p = A*(sparse(diag(1./sigx))*A');
nikcleju@62 140 opts.POSDEF = true; opts.SYM = true;
nikcleju@62 141 [dv,hcond] = linsolve(H11p, w1p, opts);
nikcleju@62 142 if (hcond < 1e-14)
nikcleju@62 143 disp('Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)');
nikcleju@62 144 xp = x;
nikcleju@62 145 return
nikcleju@62 146 end
nikcleju@62 147 dx = (w1 - w2.*sig2./sig1 - A'*dv)./sigx;
nikcleju@62 148 Adx = A*dx;
nikcleju@62 149 Atdv = A'*dv;
nikcleju@62 150 end
nikcleju@62 151
nikcleju@62 152 du = (w2 - sig2.*dx)./sig1;
nikcleju@62 153
nikcleju@62 154 dlamu1 = (lamu1./fu1).*(-dx+du) - lamu1 - (1/tau)*1./fu1;
nikcleju@62 155 dlamu2 = (lamu2./fu2).*(dx+du) - lamu2 - 1/tau*1./fu2;
nikcleju@62 156
nikcleju@62 157 % make sure that the step is feasible: keeps lamu1,lamu2 > 0, fu1,fu2 < 0
nikcleju@62 158 indp = find(dlamu1 < 0); indn = find(dlamu2 < 0);
nikcleju@62 159 s = min([1; -lamu1(indp)./dlamu1(indp); -lamu2(indn)./dlamu2(indn)]);
nikcleju@62 160 indp = find((dx-du) > 0); indn = find((-dx-du) > 0);
nikcleju@62 161 s = (0.99)*min([s; -fu1(indp)./(dx(indp)-du(indp)); -fu2(indn)./(-dx(indn)-du(indn))]);
nikcleju@62 162
nikcleju@62 163 % backtracking line search
nikcleju@62 164 suffdec = 0;
nikcleju@62 165 backiter = 0;
nikcleju@62 166 while (~suffdec)
nikcleju@62 167 xp = x + s*dx; up = u + s*du;
nikcleju@62 168 vp = v + s*dv; Atvp = Atv + s*Atdv;
nikcleju@62 169 lamu1p = lamu1 + s*dlamu1; lamu2p = lamu2 + s*dlamu2;
nikcleju@62 170 fu1p = xp - up; fu2p = -xp - up;
nikcleju@62 171 rdp = gradf0 + [lamu1p-lamu2p; -lamu1p-lamu2p] + [Atvp; zeros(N,1)];
nikcleju@62 172 rcp = [-lamu1p.*fu1p; -lamu2p.*fu2p] - (1/tau);
nikcleju@62 173 rpp = rpri + s*Adx;
nikcleju@62 174 suffdec = (norm([rdp; rcp; rpp]) <= (1-alpha*s)*resnorm);
nikcleju@62 175 s = beta*s;
nikcleju@62 176 backiter = backiter + 1;
nikcleju@62 177 if (backiter > 32)
nikcleju@62 178 disp('Stuck backtracking, returning last iterate. (See Section 4 of notes for more information.)')
nikcleju@62 179 xp = x;
nikcleju@62 180 return
nikcleju@62 181 end
nikcleju@62 182 end
nikcleju@62 183
nikcleju@62 184
nikcleju@62 185 % next iteration
nikcleju@62 186 x = xp; u = up;
nikcleju@62 187 v = vp; Atv = Atvp;
nikcleju@62 188 lamu1 = lamu1p; lamu2 = lamu2p;
nikcleju@62 189 fu1 = fu1p; fu2 = fu2p;
nikcleju@62 190
nikcleju@62 191 % surrogate duality gap
nikcleju@62 192 sdg = -(fu1'*lamu1 + fu2'*lamu2);
nikcleju@62 193 tau = mu*2*N/sdg;
nikcleju@62 194 rpri = rpp;
nikcleju@62 195 rcent = [-lamu1.*fu1; -lamu2.*fu2] - (1/tau);
nikcleju@62 196 rdual = gradf0 + [lamu1-lamu2; -lamu1-lamu2] + [Atv; zeros(N,1)];
nikcleju@62 197 resnorm = norm([rdual; rcent; rpri]);
nikcleju@62 198
nikcleju@62 199 done = (sdg < pdtol) | (pditer >= pdmaxiter);
nikcleju@62 200
nikcleju@62 201 disp(sprintf('Iteration = %d, tau = %8.3e, Primal = %8.3e, PDGap = %8.3e, Dual res = %8.3e, Primal res = %8.3e',...
nikcleju@62 202 pditer, tau, sum(u), sdg, norm(rdual), norm(rpri)));
nikcleju@62 203 if (largescale)
nikcleju@62 204 disp(sprintf(' CG Res = %8.3e, CG Iter = %d', cgres, cgiter));
nikcleju@62 205 else
nikcleju@62 206 disp(sprintf(' H11p condition number = %8.3e', hcond));
nikcleju@62 207 end
nikcleju@62 208
nikcleju@62 209 end