nikcleju@62
|
1
|
nikcleju@62
|
2 import numpy
|
nikcleju@62
|
3 import scipy.linalg
|
nikcleju@62
|
4 import math
|
nikcleju@62
|
5
|
nikcleju@63
|
6 class l1eqNotImplementedError(Exception):
|
nikcleju@62
|
7 pass
|
nikcleju@62
|
8
|
nikcleju@62
|
9 #function xp = l1eq_pd(x0, A, At, b, pdtol, pdmaxiter, cgtol, cgmaxiter)
|
nikcleju@62
|
10 def l1eq_pd(x0, A, At, b, pdtol=1e-3, pdmaxiter=50, cgtol=1e-8, cgmaxiter=200, verbose=False):
|
nikcleju@62
|
11
|
nikcleju@62
|
12 # Solve
|
nikcleju@62
|
13 # min_x ||x||_1 s.t. Ax = b
|
nikcleju@62
|
14 #
|
nikcleju@62
|
15 # Recast as linear program
|
nikcleju@62
|
16 # min_{x,u} sum(u) s.t. -u <= x <= u, Ax=b
|
nikcleju@62
|
17 # and use primal-dual interior point method
|
nikcleju@62
|
18 #
|
nikcleju@62
|
19 # Usage: xp = l1eq_pd(x0, A, At, b, pdtol, pdmaxiter, cgtol, cgmaxiter)
|
nikcleju@62
|
20 #
|
nikcleju@62
|
21 # x0 - Nx1 vector, initial point.
|
nikcleju@62
|
22 #
|
nikcleju@62
|
23 # A - Either a handle to a function that takes a N vector and returns a K
|
nikcleju@62
|
24 # vector , or a KxN matrix. If A is a function handle, the algorithm
|
nikcleju@62
|
25 # operates in "largescale" mode, solving the Newton systems via the
|
nikcleju@62
|
26 # Conjugate Gradients algorithm.
|
nikcleju@62
|
27 #
|
nikcleju@62
|
28 # At - Handle to a function that takes a K vector and returns an N vector.
|
nikcleju@62
|
29 # If A is a KxN matrix, At is ignored.
|
nikcleju@62
|
30 #
|
nikcleju@62
|
31 # b - Kx1 vector of observations.
|
nikcleju@62
|
32 #
|
nikcleju@62
|
33 # pdtol - Tolerance for primal-dual algorithm (algorithm terminates if
|
nikcleju@62
|
34 # the duality gap is less than pdtol).
|
nikcleju@62
|
35 # Default = 1e-3.
|
nikcleju@62
|
36 #
|
nikcleju@62
|
37 # pdmaxiter - Maximum number of primal-dual iterations.
|
nikcleju@62
|
38 # Default = 50.
|
nikcleju@62
|
39 #
|
nikcleju@62
|
40 # cgtol - Tolerance for Conjugate Gradients; ignored if A is a matrix.
|
nikcleju@62
|
41 # Default = 1e-8.
|
nikcleju@62
|
42 #
|
nikcleju@62
|
43 # cgmaxiter - Maximum number of iterations for Conjugate Gradients; ignored
|
nikcleju@62
|
44 # if A is a matrix.
|
nikcleju@62
|
45 # Default = 200.
|
nikcleju@62
|
46 #
|
nikcleju@62
|
47 # Written by: Justin Romberg, Caltech
|
nikcleju@62
|
48 # Email: jrom@acm.caltech.edu
|
nikcleju@62
|
49 # Created: October 2005
|
nikcleju@62
|
50
|
nikcleju@62
|
51
|
nikcleju@62
|
52 #---------------------
|
nikcleju@62
|
53 # Original Matab code:
|
nikcleju@62
|
54
|
nikcleju@62
|
55 #largescale = isa(A,'function_handle');
|
nikcleju@62
|
56 #
|
nikcleju@62
|
57 #if (nargin < 5), pdtol = 1e-3; end
|
nikcleju@62
|
58 #if (nargin < 6), pdmaxiter = 50; end
|
nikcleju@62
|
59 #if (nargin < 7), cgtol = 1e-8; end
|
nikcleju@62
|
60 #if (nargin < 8), cgmaxiter = 200; end
|
nikcleju@62
|
61 #
|
nikcleju@62
|
62 #N = length(x0);
|
nikcleju@62
|
63 #
|
nikcleju@62
|
64 #alpha = 0.01;
|
nikcleju@62
|
65 #beta = 0.5;
|
nikcleju@62
|
66 #mu = 10;
|
nikcleju@62
|
67 #
|
nikcleju@62
|
68 #gradf0 = [zeros(N,1); ones(N,1)];
|
nikcleju@62
|
69 #
|
nikcleju@62
|
70 ## starting point --- make sure that it is feasible
|
nikcleju@62
|
71 #if (largescale)
|
nikcleju@62
|
72 # if (norm(A(x0)-b)/norm(b) > cgtol)
|
nikcleju@62
|
73 # disp('Starting point infeasible; using x0 = At*inv(AAt)*y.');
|
nikcleju@62
|
74 # AAt = @(z) A(At(z));
|
nikcleju@62
|
75 # [w, cgres, cgiter] = cgsolve(AAt, b, cgtol, cgmaxiter, 0);
|
nikcleju@62
|
76 # if (cgres > 1/2)
|
nikcleju@62
|
77 # disp('A*At is ill-conditioned: cannot find starting point');
|
nikcleju@62
|
78 # xp = x0;
|
nikcleju@62
|
79 # return;
|
nikcleju@62
|
80 # end
|
nikcleju@62
|
81 # x0 = At(w);
|
nikcleju@62
|
82 # end
|
nikcleju@62
|
83 #else
|
nikcleju@62
|
84 # if (norm(A*x0-b)/norm(b) > cgtol)
|
nikcleju@62
|
85 # disp('Starting point infeasible; using x0 = At*inv(AAt)*y.');
|
nikcleju@62
|
86 # opts.POSDEF = true; opts.SYM = true;
|
nikcleju@62
|
87 # [w, hcond] = linsolve(A*A', b, opts);
|
nikcleju@62
|
88 # if (hcond < 1e-14)
|
nikcleju@62
|
89 # disp('A*At is ill-conditioned: cannot find starting point');
|
nikcleju@62
|
90 # xp = x0;
|
nikcleju@62
|
91 # return;
|
nikcleju@62
|
92 # end
|
nikcleju@62
|
93 # x0 = A'*w;
|
nikcleju@62
|
94 # end
|
nikcleju@62
|
95 #end
|
nikcleju@62
|
96 #x = x0;
|
nikcleju@62
|
97 #u = (0.95)*abs(x0) + (0.10)*max(abs(x0));
|
nikcleju@62
|
98 #
|
nikcleju@62
|
99 ## set up for the first iteration
|
nikcleju@62
|
100 #fu1 = x - u;
|
nikcleju@62
|
101 #fu2 = -x - u;
|
nikcleju@62
|
102 #lamu1 = -1./fu1;
|
nikcleju@62
|
103 #lamu2 = -1./fu2;
|
nikcleju@62
|
104 #if (largescale)
|
nikcleju@62
|
105 # v = -A(lamu1-lamu2);
|
nikcleju@62
|
106 # Atv = At(v);
|
nikcleju@62
|
107 # rpri = A(x) - b;
|
nikcleju@62
|
108 #else
|
nikcleju@62
|
109 # v = -A*(lamu1-lamu2);
|
nikcleju@62
|
110 # Atv = A'*v;
|
nikcleju@62
|
111 # rpri = A*x - b;
|
nikcleju@62
|
112 #end
|
nikcleju@62
|
113 #
|
nikcleju@62
|
114 #sdg = -(fu1'*lamu1 + fu2'*lamu2);
|
nikcleju@62
|
115 #tau = mu*2*N/sdg;
|
nikcleju@62
|
116 #
|
nikcleju@62
|
117 #rcent = [-lamu1.*fu1; -lamu2.*fu2] - (1/tau);
|
nikcleju@62
|
118 #rdual = gradf0 + [lamu1-lamu2; -lamu1-lamu2] + [Atv; zeros(N,1)];
|
nikcleju@62
|
119 #resnorm = norm([rdual; rcent; rpri]);
|
nikcleju@62
|
120 #
|
nikcleju@62
|
121 #pditer = 0;
|
nikcleju@62
|
122 #done = (sdg < pdtol) | (pditer >= pdmaxiter);
|
nikcleju@62
|
123 #while (~done)
|
nikcleju@62
|
124 #
|
nikcleju@62
|
125 # pditer = pditer + 1;
|
nikcleju@62
|
126 #
|
nikcleju@62
|
127 # w1 = -1/tau*(-1./fu1 + 1./fu2) - Atv;
|
nikcleju@62
|
128 # w2 = -1 - 1/tau*(1./fu1 + 1./fu2);
|
nikcleju@62
|
129 # w3 = -rpri;
|
nikcleju@62
|
130 #
|
nikcleju@62
|
131 # sig1 = -lamu1./fu1 - lamu2./fu2;
|
nikcleju@62
|
132 # sig2 = lamu1./fu1 - lamu2./fu2;
|
nikcleju@62
|
133 # sigx = sig1 - sig2.^2./sig1;
|
nikcleju@62
|
134 #
|
nikcleju@62
|
135 # if (largescale)
|
nikcleju@62
|
136 # w1p = w3 - A(w1./sigx - w2.*sig2./(sigx.*sig1));
|
nikcleju@62
|
137 # h11pfun = @(z) -A(1./sigx.*At(z));
|
nikcleju@62
|
138 # [dv, cgres, cgiter] = cgsolve(h11pfun, w1p, cgtol, cgmaxiter, 0);
|
nikcleju@62
|
139 # if (cgres > 1/2)
|
nikcleju@62
|
140 # disp('Cannot solve system. Returning previous iterate. (See Section 4 of notes for more information.)');
|
nikcleju@62
|
141 # xp = x;
|
nikcleju@62
|
142 # return
|
nikcleju@62
|
143 # end
|
nikcleju@62
|
144 # dx = (w1 - w2.*sig2./sig1 - At(dv))./sigx;
|
nikcleju@62
|
145 # Adx = A(dx);
|
nikcleju@62
|
146 # Atdv = At(dv);
|
nikcleju@62
|
147 # else
|
nikcleju@62
|
148 # w1p = -(w3 - A*(w1./sigx - w2.*sig2./(sigx.*sig1)));
|
nikcleju@62
|
149 # H11p = A*(sparse(diag(1./sigx))*A');
|
nikcleju@62
|
150 # opts.POSDEF = true; opts.SYM = true;
|
nikcleju@62
|
151 # [dv,hcond] = linsolve(H11p, w1p, opts);
|
nikcleju@62
|
152 # if (hcond < 1e-14)
|
nikcleju@62
|
153 # disp('Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)');
|
nikcleju@62
|
154 # xp = x;
|
nikcleju@62
|
155 # return
|
nikcleju@62
|
156 # end
|
nikcleju@62
|
157 # dx = (w1 - w2.*sig2./sig1 - A'*dv)./sigx;
|
nikcleju@62
|
158 # Adx = A*dx;
|
nikcleju@62
|
159 # Atdv = A'*dv;
|
nikcleju@62
|
160 # end
|
nikcleju@62
|
161 #
|
nikcleju@62
|
162 # du = (w2 - sig2.*dx)./sig1;
|
nikcleju@62
|
163 #
|
nikcleju@62
|
164 # dlamu1 = (lamu1./fu1).*(-dx+du) - lamu1 - (1/tau)*1./fu1;
|
nikcleju@62
|
165 # dlamu2 = (lamu2./fu2).*(dx+du) - lamu2 - 1/tau*1./fu2;
|
nikcleju@62
|
166 #
|
nikcleju@62
|
167 # # make sure that the step is feasible: keeps lamu1,lamu2 > 0, fu1,fu2 < 0
|
nikcleju@62
|
168 # indp = find(dlamu1 < 0); indn = find(dlamu2 < 0);
|
nikcleju@62
|
169 # s = min([1; -lamu1(indp)./dlamu1(indp); -lamu2(indn)./dlamu2(indn)]);
|
nikcleju@62
|
170 # indp = find((dx-du) > 0); indn = find((-dx-du) > 0);
|
nikcleju@62
|
171 # s = (0.99)*min([s; -fu1(indp)./(dx(indp)-du(indp)); -fu2(indn)./(-dx(indn)-du(indn))]);
|
nikcleju@62
|
172 #
|
nikcleju@62
|
173 # # backtracking line search
|
nikcleju@62
|
174 # suffdec = 0;
|
nikcleju@62
|
175 # backiter = 0;
|
nikcleju@62
|
176 # while (~suffdec)
|
nikcleju@62
|
177 # xp = x + s*dx; up = u + s*du;
|
nikcleju@62
|
178 # vp = v + s*dv; Atvp = Atv + s*Atdv;
|
nikcleju@62
|
179 # lamu1p = lamu1 + s*dlamu1; lamu2p = lamu2 + s*dlamu2;
|
nikcleju@62
|
180 # fu1p = xp - up; fu2p = -xp - up;
|
nikcleju@62
|
181 # rdp = gradf0 + [lamu1p-lamu2p; -lamu1p-lamu2p] + [Atvp; zeros(N,1)];
|
nikcleju@62
|
182 # rcp = [-lamu1p.*fu1p; -lamu2p.*fu2p] - (1/tau);
|
nikcleju@62
|
183 # rpp = rpri + s*Adx;
|
nikcleju@62
|
184 # suffdec = (norm([rdp; rcp; rpp]) <= (1-alpha*s)*resnorm);
|
nikcleju@62
|
185 # s = beta*s;
|
nikcleju@62
|
186 # backiter = backiter + 1;
|
nikcleju@62
|
187 # if (backiter > 32)
|
nikcleju@62
|
188 # disp('Stuck backtracking, returning last iterate. (See Section 4 of notes for more information.)')
|
nikcleju@62
|
189 # xp = x;
|
nikcleju@62
|
190 # return
|
nikcleju@62
|
191 # end
|
nikcleju@62
|
192 # end
|
nikcleju@62
|
193 #
|
nikcleju@62
|
194 #
|
nikcleju@62
|
195 # # next iteration
|
nikcleju@62
|
196 # x = xp; u = up;
|
nikcleju@62
|
197 # v = vp; Atv = Atvp;
|
nikcleju@62
|
198 # lamu1 = lamu1p; lamu2 = lamu2p;
|
nikcleju@62
|
199 # fu1 = fu1p; fu2 = fu2p;
|
nikcleju@62
|
200 #
|
nikcleju@62
|
201 # # surrogate duality gap
|
nikcleju@62
|
202 # sdg = -(fu1'*lamu1 + fu2'*lamu2);
|
nikcleju@62
|
203 # tau = mu*2*N/sdg;
|
nikcleju@62
|
204 # rpri = rpp;
|
nikcleju@62
|
205 # rcent = [-lamu1.*fu1; -lamu2.*fu2] - (1/tau);
|
nikcleju@62
|
206 # rdual = gradf0 + [lamu1-lamu2; -lamu1-lamu2] + [Atv; zeros(N,1)];
|
nikcleju@62
|
207 # resnorm = norm([rdual; rcent; rpri]);
|
nikcleju@62
|
208 #
|
nikcleju@62
|
209 # done = (sdg < pdtol) | (pditer >= pdmaxiter);
|
nikcleju@62
|
210 #
|
nikcleju@62
|
211 # disp(sprintf('Iteration = #d, tau = #8.3e, Primal = #8.3e, PDGap = #8.3e, Dual res = #8.3e, Primal res = #8.3e',...
|
nikcleju@62
|
212 # pditer, tau, sum(u), sdg, norm(rdual), norm(rpri)));
|
nikcleju@62
|
213 # if (largescale)
|
nikcleju@62
|
214 # disp(sprintf(' CG Res = #8.3e, CG Iter = #d', cgres, cgiter));
|
nikcleju@62
|
215 # else
|
nikcleju@62
|
216 # disp(sprintf(' H11p condition number = #8.3e', hcond));
|
nikcleju@62
|
217 # end
|
nikcleju@62
|
218 #
|
nikcleju@62
|
219 #end
|
nikcleju@62
|
220
|
nikcleju@62
|
221 # End of original Matab code
|
nikcleju@62
|
222 #----------------------------
|
nikcleju@62
|
223
|
nikcleju@65
|
224 # Nic: check if b is 0; if so, return 0
|
nikcleju@65
|
225 # Otherwise it will break later
|
nikcleju@65
|
226 if numpy.linalg.norm(b) < 1e-16:
|
nikcleju@65
|
227 return numpy.zeros_like(x0)
|
nikcleju@65
|
228
|
nikcleju@62
|
229 #largescale = isa(A,'function_handle');
|
nikcleju@62
|
230 if hasattr(A, '__call__'):
|
nikcleju@62
|
231 largescale = True
|
nikcleju@62
|
232 else:
|
nikcleju@62
|
233 largescale = False
|
nikcleju@62
|
234
|
nikcleju@62
|
235 #N = length(x0);
|
nikcleju@62
|
236 N = x0.size
|
nikcleju@62
|
237
|
nikcleju@62
|
238 alpha = 0.01
|
nikcleju@62
|
239 beta = 0.5
|
nikcleju@62
|
240 mu = 10
|
nikcleju@62
|
241
|
nikcleju@62
|
242 #gradf0 = [zeros(N,1); ones(N,1)];
|
nikcleju@62
|
243 gradf0 = numpy.hstack((numpy.zeros(N), numpy.ones(N)))
|
nikcleju@62
|
244
|
nikcleju@62
|
245 # starting point --- make sure that it is feasible
|
nikcleju@62
|
246 #if (largescale)
|
nikcleju@62
|
247 if largescale:
|
nikcleju@63
|
248 raise l1eqNotImplementedError('Largescale not implemented yet!')
|
nikcleju@62
|
249 else:
|
nikcleju@62
|
250 #if (norm(A*x0-b)/norm(b) > cgtol)
|
nikcleju@62
|
251 if numpy.linalg.norm(numpy.dot(A,x0)-b) / numpy.linalg.norm(b) > cgtol:
|
nikcleju@62
|
252 #disp('Starting point infeasible; using x0 = At*inv(AAt)*y.');
|
nikcleju@62
|
253 if verbose:
|
nikcleju@62
|
254 print 'Starting point infeasible; using x0 = At*inv(AAt)*y.'
|
nikcleju@62
|
255 #opts.POSDEF = true; opts.SYM = true;
|
nikcleju@62
|
256 #[w, hcond] = linsolve(A*A', b, opts);
|
nikcleju@62
|
257 #if (hcond < 1e-14)
|
nikcleju@62
|
258 # disp('A*At is ill-conditioned: cannot find starting point');
|
nikcleju@62
|
259 # xp = x0;
|
nikcleju@62
|
260 # return;
|
nikcleju@62
|
261 #end
|
nikcleju@62
|
262 #x0 = A'*w;
|
nikcleju@62
|
263 try:
|
nikcleju@62
|
264 w = scipy.linalg.solve(numpy.dot(A,A.T), b, sym_pos=True)
|
nikcleju@65
|
265 hcond = 1.0/numpy.linalg.cond(numpy.dot(A,A.T))
|
nikcleju@62
|
266 except scipy.linalg.LinAlgError:
|
nikcleju@62
|
267 if verbose:
|
nikcleju@62
|
268 print 'A*At is ill-conditioned: cannot find starting point'
|
nikcleju@62
|
269 xp = x0.copy()
|
nikcleju@62
|
270 return xp
|
nikcleju@62
|
271 if hcond < 1e-14:
|
nikcleju@62
|
272 if verbose:
|
nikcleju@62
|
273 print 'A*At is ill-conditioned: cannot find starting point'
|
nikcleju@62
|
274 xp = x0.copy()
|
nikcleju@62
|
275 return xp
|
nikcleju@62
|
276 x0 = numpy.dot(A.T, w)
|
nikcleju@62
|
277 #end
|
nikcleju@62
|
278 #end
|
nikcleju@62
|
279 x = x0.copy()
|
nikcleju@62
|
280 #u = (0.95)*abs(x0) + (0.10)*max(abs(x0));
|
nikcleju@62
|
281 u = (0.95)*numpy.abs(x0) + (0.10)*numpy.abs(x0).max()
|
nikcleju@62
|
282
|
nikcleju@62
|
283 # set up for the first iteration
|
nikcleju@62
|
284 fu1 = x - u
|
nikcleju@62
|
285 fu2 = -x - u
|
nikcleju@62
|
286 lamu1 = -1/fu1
|
nikcleju@62
|
287 lamu2 = -1/fu2
|
nikcleju@62
|
288 if (largescale):
|
nikcleju@62
|
289 #v = -A(lamu1-lamu2);
|
nikcleju@62
|
290 #Atv = At(v);
|
nikcleju@62
|
291 #rpri = A(x) - b;
|
nikcleju@63
|
292 raise l1eqNotImplementedError('Largescale not implemented yet!')
|
nikcleju@62
|
293 else:
|
nikcleju@62
|
294 #v = -A*(lamu1-lamu2);
|
nikcleju@62
|
295 #Atv = A'*v;
|
nikcleju@62
|
296 #rpri = A*x - b;
|
nikcleju@62
|
297 v = numpy.dot(-A, lamu1-lamu2)
|
nikcleju@62
|
298 Atv = numpy.dot(A.T, v)
|
nikcleju@62
|
299 rpri = numpy.dot(A,x) - b
|
nikcleju@62
|
300 #end
|
nikcleju@62
|
301
|
nikcleju@62
|
302 #sdg = -(fu1'*lamu1 + fu2'*lamu2);
|
nikcleju@62
|
303 sdg = -(numpy.dot(fu1,lamu1) + numpy.dot(fu2,lamu2))
|
nikcleju@62
|
304 tau = mu*2*N/sdg
|
nikcleju@62
|
305
|
nikcleju@62
|
306 #rcent = [-lamu1.*fu1; -lamu2.*fu2] - (1/tau);
|
nikcleju@62
|
307 rcent = numpy.hstack((-numpy.dot(lamu1,fu1), -numpy.dot(lamu2,fu2))) - (1/tau)
|
nikcleju@62
|
308 #rdual = gradf0 + [lamu1-lamu2; -lamu1-lamu2] + [Atv; zeros(N,1)];
|
nikcleju@62
|
309 rdual = gradf0 + numpy.hstack((lamu1-lamu2, -lamu1-lamu2)) + numpy.hstack((Atv, numpy.zeros(N)))
|
nikcleju@62
|
310 #resnorm = norm([rdual; rcent; rpri]);
|
nikcleju@62
|
311 resnorm = numpy.linalg.norm(numpy.hstack((rdual, rcent, rpri)))
|
nikcleju@62
|
312
|
nikcleju@62
|
313 pditer = 0
|
nikcleju@62
|
314 #done = (sdg < pdtol) | (pditer >= pdmaxiter);
|
nikcleju@62
|
315 done = (sdg < pdtol) or (pditer >= pdmaxiter)
|
nikcleju@62
|
316 #while (~done)
|
nikcleju@62
|
317 while not done:
|
nikcleju@62
|
318
|
nikcleju@62
|
319 pditer = pditer + 1
|
nikcleju@62
|
320
|
nikcleju@62
|
321 #w1 = -1/tau*(-1./fu1 + 1./fu2) - Atv;
|
nikcleju@62
|
322 w1 = -1/tau*(-1/fu1 + 1/fu2) - Atv
|
nikcleju@62
|
323 w2 = -1 - 1/tau*(1/fu1 + 1/fu2)
|
nikcleju@62
|
324 w3 = -rpri
|
nikcleju@62
|
325
|
nikcleju@62
|
326 #sig1 = -lamu1./fu1 - lamu2./fu2;
|
nikcleju@62
|
327 sig1 = -lamu1/fu1 - lamu2/fu2
|
nikcleju@62
|
328 sig2 = lamu1/fu1 - lamu2/fu2
|
nikcleju@62
|
329 #sigx = sig1 - sig2.^2./sig1;
|
nikcleju@62
|
330 sigx = sig1 - sig2**2/sig1
|
nikcleju@62
|
331
|
nikcleju@62
|
332 if largescale:
|
nikcleju@62
|
333 #w1p = w3 - A(w1./sigx - w2.*sig2./(sigx.*sig1));
|
nikcleju@62
|
334 #h11pfun = @(z) -A(1./sigx.*At(z));
|
nikcleju@62
|
335 #[dv, cgres, cgiter] = cgsolve(h11pfun, w1p, cgtol, cgmaxiter, 0);
|
nikcleju@62
|
336 #if (cgres > 1/2)
|
nikcleju@62
|
337 # disp('Cannot solve system. Returning previous iterate. (See Section 4 of notes for more information.)');
|
nikcleju@62
|
338 # xp = x;
|
nikcleju@62
|
339 # return
|
nikcleju@62
|
340 #end
|
nikcleju@62
|
341 #dx = (w1 - w2.*sig2./sig1 - At(dv))./sigx;
|
nikcleju@62
|
342 #Adx = A(dx);
|
nikcleju@62
|
343 #Atdv = At(dv);
|
nikcleju@63
|
344 raise l1eqNotImplementedError('Largescale not implemented yet!')
|
nikcleju@62
|
345 else:
|
nikcleju@62
|
346 #w1p = -(w3 - A*(w1./sigx - w2.*sig2./(sigx.*sig1)));
|
nikcleju@62
|
347 w1p = -(w3 - numpy.dot(A,(w1/sigx - w2*sig2/(sigx*sig1))))
|
nikcleju@62
|
348 #H11p = A*(sparse(diag(1./sigx))*A');
|
nikcleju@62
|
349 H11p = numpy.dot(A, numpy.dot(numpy.diag(1/sigx),A.T))
|
nikcleju@62
|
350 #opts.POSDEF = true; opts.SYM = true;
|
nikcleju@62
|
351 #[dv,hcond] = linsolve(H11p, w1p, opts);
|
nikcleju@62
|
352 try:
|
nikcleju@62
|
353 dv = scipy.linalg.solve(H11p, w1p, sym_pos=True)
|
nikcleju@62
|
354 hcond = 1.0/numpy.linalg.cond(H11p)
|
nikcleju@62
|
355 except scipy.linalg.LinAlgError:
|
nikcleju@62
|
356 if verbose:
|
nikcleju@62
|
357 print 'Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)'
|
nikcleju@62
|
358 xp = x.copy()
|
nikcleju@62
|
359 return xp
|
nikcleju@62
|
360 if hcond < 1e-14:
|
nikcleju@62
|
361 if verbose:
|
nikcleju@62
|
362 print 'Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)'
|
nikcleju@62
|
363 xp = x.copy()
|
nikcleju@62
|
364 return xp
|
nikcleju@62
|
365 #if (hcond < 1e-14)
|
nikcleju@62
|
366 # disp('Matrix ill-conditioned. Returning previous iterate. (See Section 4 of notes for more information.)');
|
nikcleju@62
|
367 # xp = x;
|
nikcleju@62
|
368 # return
|
nikcleju@62
|
369 #end
|
nikcleju@62
|
370
|
nikcleju@62
|
371 #dx = (w1 - w2.*sig2./sig1 - A'*dv)./sigx;
|
nikcleju@62
|
372 dx = (w1 - w2*sig2/sig1 - numpy.dot(A.T,dv))/sigx
|
nikcleju@62
|
373 #Adx = A*dx;
|
nikcleju@62
|
374 Adx = numpy.dot(A,dx)
|
nikcleju@62
|
375 #Atdv = A'*dv;
|
nikcleju@62
|
376 Atdv = numpy.dot(A.T,dv)
|
nikcleju@62
|
377 #end
|
nikcleju@62
|
378
|
nikcleju@62
|
379 #du = (w2 - sig2.*dx)./sig1;
|
nikcleju@62
|
380 du = (w2 - sig2*dx)/sig1
|
nikcleju@62
|
381
|
nikcleju@62
|
382 #dlamu1 = (lamu1./fu1).*(-dx+du) - lamu1 - (1/tau)*1./fu1;
|
nikcleju@62
|
383 dlamu1 = (lamu1/fu1)*(-dx+du) - lamu1 - (1/tau)*1/fu1
|
nikcleju@62
|
384 dlamu2 = (lamu2/fu2)*(dx+du) - lamu2 - 1/tau*1/fu2
|
nikcleju@62
|
385
|
nikcleju@62
|
386 # make sure that the step is feasible: keeps lamu1,lamu2 > 0, fu1,fu2 < 0
|
nikcleju@62
|
387 #indp = find(dlamu1 < 0); indn = find(dlamu2 < 0);
|
nikcleju@62
|
388 indp = numpy.nonzero(dlamu1 < 0)
|
nikcleju@62
|
389 indn = numpy.nonzero(dlamu2 < 0)
|
nikcleju@62
|
390 #s = min([1; -lamu1(indp)./dlamu1(indp); -lamu2(indn)./dlamu2(indn)]);
|
nikcleju@62
|
391 s = numpy.min(numpy.hstack((numpy.array([1]), -lamu1[indp]/dlamu1[indp], -lamu2[indn]/dlamu2[indn])))
|
nikcleju@62
|
392 #indp = find((dx-du) > 0); indn = find((-dx-du) > 0);
|
nikcleju@62
|
393 indp = numpy.nonzero((dx-du) > 0)
|
nikcleju@62
|
394 indn = numpy.nonzero((-dx-du) > 0)
|
nikcleju@62
|
395 #s = (0.99)*min([s; -fu1(indp)./(dx(indp)-du(indp)); -fu2(indn)./(-dx(indn)-du(indn))]);
|
nikcleju@62
|
396 s = (0.99)*numpy.min(numpy.hstack((numpy.array([s]), -fu1[indp]/(dx[indp]-du[indp]), -fu2[indn]/(-dx[indn]-du[indn]))))
|
nikcleju@62
|
397
|
nikcleju@62
|
398 # backtracking line search
|
nikcleju@62
|
399 suffdec = 0
|
nikcleju@62
|
400 backiter = 0
|
nikcleju@62
|
401 #while (~suffdec)
|
nikcleju@62
|
402 while not suffdec:
|
nikcleju@62
|
403 #xp = x + s*dx; up = u + s*du;
|
nikcleju@62
|
404 xp = x + s*dx
|
nikcleju@62
|
405 up = u + s*du
|
nikcleju@62
|
406 #vp = v + s*dv; Atvp = Atv + s*Atdv;
|
nikcleju@62
|
407 vp = v + s*dv
|
nikcleju@62
|
408 Atvp = Atv + s*Atdv
|
nikcleju@62
|
409 #lamu1p = lamu1 + s*dlamu1; lamu2p = lamu2 + s*dlamu2;
|
nikcleju@62
|
410 lamu1p = lamu1 + s*dlamu1
|
nikcleju@62
|
411 lamu2p = lamu2 + s*dlamu2
|
nikcleju@62
|
412 #fu1p = xp - up; fu2p = -xp - up;
|
nikcleju@62
|
413 fu1p = xp - up
|
nikcleju@62
|
414 fu2p = -xp - up
|
nikcleju@62
|
415 #rdp = gradf0 + [lamu1p-lamu2p; -lamu1p-lamu2p] + [Atvp; zeros(N,1)];
|
nikcleju@62
|
416 rdp = gradf0 + numpy.hstack((lamu1p-lamu2p, -lamu1p-lamu2p)) + numpy.hstack((Atvp, numpy.zeros(N)))
|
nikcleju@62
|
417 #rcp = [-lamu1p.*fu1p; -lamu2p.*fu2p] - (1/tau);
|
nikcleju@62
|
418 rcp = numpy.hstack((-lamu1p*fu1p, -lamu2p*fu2p)) - (1/tau)
|
nikcleju@62
|
419 #rpp = rpri + s*Adx;
|
nikcleju@62
|
420 rpp = rpri + s*Adx
|
nikcleju@62
|
421 #suffdec = (norm([rdp; rcp; rpp]) <= (1-alpha*s)*resnorm);
|
nikcleju@62
|
422 suffdec = (numpy.linalg.norm(numpy.hstack((rdp, rcp, rpp))) <= (1-alpha*s)*resnorm)
|
nikcleju@62
|
423 s = beta*s
|
nikcleju@62
|
424 backiter = backiter + 1
|
nikcleju@62
|
425 if (backiter > 32):
|
nikcleju@62
|
426 if verbose:
|
nikcleju@62
|
427 print 'Stuck backtracking, returning last iterate. (See Section 4 of notes for more information.)'
|
nikcleju@62
|
428 xp = x.copy()
|
nikcleju@62
|
429 return xp
|
nikcleju@62
|
430 #end
|
nikcleju@62
|
431 #end
|
nikcleju@62
|
432
|
nikcleju@62
|
433
|
nikcleju@62
|
434 # next iteration
|
nikcleju@62
|
435 #x = xp; u = up;
|
nikcleju@62
|
436 x = xp.copy()
|
nikcleju@62
|
437 u = up.copy()
|
nikcleju@62
|
438 #v = vp; Atv = Atvp;
|
nikcleju@62
|
439 v = vp.copy()
|
nikcleju@62
|
440 Atv = Atvp.copy()
|
nikcleju@62
|
441 #lamu1 = lamu1p; lamu2 = lamu2p;
|
nikcleju@62
|
442 lamu1 = lamu1p.copy()
|
nikcleju@62
|
443 lamu2 = lamu2p.copy()
|
nikcleju@62
|
444 #fu1 = fu1p; fu2 = fu2p;
|
nikcleju@62
|
445 fu1 = fu1p.copy()
|
nikcleju@62
|
446 fu2 = fu2p.copy()
|
nikcleju@62
|
447
|
nikcleju@62
|
448 # surrogate duality gap
|
nikcleju@62
|
449 #sdg = -(fu1'*lamu1 + fu2'*lamu2);
|
nikcleju@62
|
450 sdg = -(numpy.dot(fu1,lamu1) + numpy.dot(fu2,lamu2))
|
nikcleju@62
|
451 tau = mu*2*N/sdg
|
nikcleju@62
|
452 rpri = rpp.copy()
|
nikcleju@62
|
453 #rcent = [-lamu1.*fu1; -lamu2.*fu2] - (1/tau);
|
nikcleju@62
|
454 rcent = numpy.hstack((-lamu1*fu1, -lamu2*fu2)) - (1/tau)
|
nikcleju@62
|
455 #rdual = gradf0 + [lamu1-lamu2; -lamu1-lamu2] + [Atv; zeros(N,1)];
|
nikcleju@62
|
456 rdual = gradf0 + numpy.hstack((lamu1-lamu2, -lamu1-lamu2)) + numpy.hstack((Atv, numpy.zeros(N)))
|
nikcleju@62
|
457 #resnorm = norm([rdual; rcent; rpri]);
|
nikcleju@62
|
458 resnorm = numpy.linalg.norm(numpy.hstack((rdual, rcent, rpri)))
|
nikcleju@62
|
459
|
nikcleju@62
|
460 #done = (sdg < pdtol) | (pditer >= pdmaxiter);
|
nikcleju@62
|
461 done = (sdg < pdtol) or (pditer >= pdmaxiter)
|
nikcleju@62
|
462
|
nikcleju@62
|
463 if verbose:
|
nikcleju@62
|
464 print 'Iteration =',pditer,', tau =',tau,', Primal =',numpy.sum(u),', PDGap =',sdg,', Dual res =',numpy.linalg.norm(rdual),', Primal res =',numpy.linalg.norm(rpri)
|
nikcleju@62
|
465 if largescale:
|
nikcleju@62
|
466 #disp(sprintf(' CG Res = #8.3e, CG Iter = #d', cgres, cgiter));
|
nikcleju@63
|
467 raise l1eqNotImplementedError('Largescale not implemented yet!')
|
nikcleju@62
|
468 else:
|
nikcleju@62
|
469 #disp(sprintf(' H11p condition number = #8.3e', hcond));
|
nikcleju@62
|
470 if verbose:
|
nikcleju@62
|
471 print ' H11p condition number =',hcond
|
nikcleju@62
|
472 #end
|
nikcleju@62
|
473
|
nikcleju@62
|
474 #end
|
nikcleju@62
|
475
|
nikcleju@62
|
476 return xp |