Mercurial > hg > smallbox
comparison toolboxes/alps/ALPS/l_ALPS.m @ 154:0de08f68256b ivand_dev
ALPS toolbox - Algebraic Pursuit added to smallbox
author | Ivan Damnjanovic lnx <ivan.damnjanovic@eecs.qmul.ac.uk> |
---|---|
date | Fri, 12 Aug 2011 11:17:47 +0100 |
parents | |
children |
comparison
equal
deleted
inserted
replaced
153:af307f247ac7 | 154:0de08f68256b |
---|---|
1 function [x_hat, numiter, x_path] = l_ALPS(y, Phi, K, params) | |
2 % ========================================================================= | |
3 % l-ALPS(#) algorithm - Beta Version | |
4 % ========================================================================= | |
5 % Algebraic Pursuit (ALPS) algorithm with l-memory acceleration. | |
6 % | |
7 % Detailed discussion on the algorithm can be found in section III in | |
8 % [1] "On Accelerated Hard Thresholding Methods for Sparse Approximation", written | |
9 % by Volkan Cevher, Technical Report, 2011. | |
10 % ========================================================================= | |
11 % INPUT ARGUMENTS: | |
12 % y M x 1 undersampled measurement vector. | |
13 % Phi M x N regression matrix. | |
14 % K Sparsity of underlying vector x* or desired | |
15 % sparsity of solution. | |
16 % params Structure of parameters. These are: | |
17 % | |
18 % tol,... Early stopping tolerance. Default value: tol = | |
19 % 1-e5 | |
20 % ALPSiters,... Maximum number of algorithm iterations. Default | |
21 % value: 300. | |
22 % solveNewtonb,... If solveNewtonb == 1: Corresponds to solving a | |
23 % Newton system restricted to a sparse support. | |
24 % It is implemented via conjugate gradients. | |
25 % If solveNewtonb == 0: Step size selection as described | |
26 % in eqs. (12) and (13) in [1]. | |
27 % Default value: solveNewtonb = 0. | |
28 % gradientDescentx,... If gradientDescentx == 1: single gradient | |
29 % update of x_{i+1} restricted ot its support with | |
30 % line search. Default value: gradientDescentx = | |
31 % 1. | |
32 % solveNewtonx,... If solveNewtonx == 1: Akin to Hard Thresholding Pursuit | |
33 % (c.f. Simon Foucart, "Hard Thresholding Pursuit," | |
34 % preprint, 2010). Default vale: solveNewtonx = 0 | |
35 % tau,... Variable that controls the momentum in | |
36 % non-memoryless case. Ignored in memoryless | |
37 % case. Default value: tau = 1/2. | |
38 % Special cases: | |
39 % - tau = -1: momentum step size is automatically | |
40 % optimized in every step. | |
41 % - tau as a function handle: user defined | |
42 % behavior of tau momentum term. | |
43 % mu,... Variable that controls the step size selection. | |
44 % When mu = 0, step size is computed adaptively | |
45 % per iteration. Default value: mu = 0. | |
46 % cg_maxiter,... Maximum iterations for Conjugate-Gradients method. | |
47 % cg_tol Tolerance variable for Conjugate-Gradients method. | |
48 % ========================================================================= | |
49 % OUTPUT ARGUMENTS: | |
50 % x_hat N x 1 recovered K-sparse vector. | |
51 % numiter Number of iterations executed. | |
52 % x_path Keeps a series of computed N x 1 K-sparse vectors | |
53 % until the end of the iterative process. | |
54 % ========================================================================= | |
55 % 01/04/2011, by Anastasios Kyrillidis. anastasios.kyrillidis@epfl.ch, EPFL. | |
56 % ========================================================================= | |
57 % cgsolve.m is written by Justin Romberg, Caltech, Oct. 2005. | |
58 % Email: jrom@acm.caltech.edu | |
59 % ========================================================================= | |
60 % This work was supported in part by the European Commission under Grant | |
61 % MIRG-268398 and DARPA KeCoM program #11-DARPA-1055. VC also would like | |
62 % to acknowledge Rice University for his Faculty Fellowship. | |
63 % ========================================================================= | |
64 | |
65 [M, N] = size(Phi); | |
66 | |
67 %% Initialize transpose of measurement matrix | |
68 | |
69 Phi_t = Phi'; | |
70 | |
71 %% Initialize to zero vector | |
72 y_cur = zeros(N,1); | |
73 Y_i = []; | |
74 | |
75 x_path = zeros(N, params.ALPSiters); | |
76 | |
77 %% CG params | |
78 if (params.solveNewtonx == 1 || params.solveNewtonb == 1) | |
79 cg_verbose = 0; | |
80 cg_A = Phi_t*Phi; | |
81 cg_b = Phi_t*y; | |
82 end; | |
83 | |
84 %% Determine momentum step size selection strategy | |
85 optimizeTau = 0; | |
86 function_tau = 0; | |
87 | |
88 if (isa(params.tau,'float')) | |
89 if (params.tau == -1) | |
90 optimizeTau = 1; | |
91 end; | |
92 elseif (isa(params.tau, 'function_handle')) | |
93 function_tau = 1; | |
94 end; | |
95 | |
96 %% Determine step size selection strategy | |
97 function_mu = 0; | |
98 adaptive_mu = 0; | |
99 | |
100 if (isa(params.mu,'float')) | |
101 function_mu = 0; | |
102 if (params.mu == 0) | |
103 adaptive_mu = 1; | |
104 else | |
105 adaptive_mu = 0; | |
106 end; | |
107 elseif (isa(params.mu,'function_handle')) | |
108 function_mu = 1; | |
109 end; | |
110 | |
111 %% Help variables | |
112 complementary_Yi = ones(N,1); | |
113 x = zeros(N,params.memory + 1); | |
114 Phi_x = zeros(M,params.memory + 1); | |
115 tau_candidates = zeros(params.memory,1); | |
116 y_candidates = zeros(N,params.memory); | |
117 residual_energy = zeros(params.memory,1); | |
118 | |
119 i = 1; | |
120 %% l-ALPS(#) | |
121 while (i <= params.ALPSiters) | |
122 x_path(:,i) = x(:,end); | |
123 | |
124 for k = 1:params.memory | |
125 x(:,k) = x(:,k+1); | |
126 Phi_x(:,k) = Phi_x(:,k+1); | |
127 end; | |
128 | |
129 % Compute the residual | |
130 if (i == 1) | |
131 res = y; | |
132 else | |
133 Phi_y_cur = Phi(:,Y_i)*y_cur(Y_i); | |
134 res = y - Phi_y_cur; | |
135 end; | |
136 | |
137 % Compute the derivative | |
138 der = Phi_t*res; | |
139 | |
140 % Determine S_i via eq. (11) (change of variable from x_i to y_i) | |
141 complementary_Yi(Y_i) = 0; | |
142 [tmpArg, ind_der] = sort(abs(der).*complementary_Yi, 'descend'); | |
143 complementary_Yi(Y_i) = 1; | |
144 S_i = [Y_i; ind_der(1:K)]; | |
145 | |
146 ider = der(S_i); | |
147 if (params.solveNewtonb == 1) | |
148 % Compute least squares solution of the system A*y = (A*A)x using CG | |
149 if (params.useCG == 1) | |
150 [b, tmpArg, tmpArg] = cgsolve(cg_A(S_i, S_i), cg_b(S_i), params.cg_tol, params.cg_maxiter, cg_verbose); | |
151 else | |
152 b = cg_A(S_i,S_i)\cg_b(S_i); | |
153 end; | |
154 else | |
155 % Step size selection via eq. (12) and eq. (13) (change of variable from x_i to y_i) | |
156 if (adaptive_mu) | |
157 Pder = Phi(:,S_i)*ider; | |
158 mu_bar = ider'*ider/(Pder'*Pder); | |
159 b = y_cur(S_i) + mu_bar*ider; | |
160 elseif (function_mu) | |
161 b = y_cur(S_i) + params.mu(i)*ider; | |
162 else b = y_cur(S_i) + params.mu*ider; | |
163 end; | |
164 end; | |
165 | |
166 % Hard-threshold b and compute X_{i+1} | |
167 [tmpArg, ind_b] = sort(abs(b), 'descend'); | |
168 x(:,end) = zeros(N,1); | |
169 x(S_i(ind_b(1:K)),end) = b(ind_b(1:K)); | |
170 X_i = S_i(ind_b(1:K)); | |
171 | |
172 if (params.gradientDescentx == 1) | |
173 % Calculate gradient of estimated current x vector | |
174 Phi_x_cur = Phi(:,X_i)*x(X_i,end); | |
175 res = y - Phi_x_cur; | |
176 der = Phi_t*res; | |
177 ider = der(X_i); | |
178 | |
179 if (adaptive_mu) | |
180 Pder = Phi(:,X_i)*ider; | |
181 mu_bar = ider'*ider/(Pder'*Pder); | |
182 x(X_i,end) = x(X_i,end) + mu_bar*ider; | |
183 elseif (function_mu) | |
184 x(X_i,end) = x(X_i,end) + params.mu(i)*ider; | |
185 else x(X_i,end) = x(X_i,end) + params.mu*ider; | |
186 end; | |
187 elseif (params.solveNewtonx == 1) | |
188 % Similar to HTP | |
189 if (params.useCG == 1) | |
190 [v, tmpArg, tmpArg] = cgsolve(cg_A(X_i, X_i), cg_b(X_i), params.cg_tol, params.cg_maxiter, cg_verbose); | |
191 else | |
192 v = cg_A(X_i, X_i)\cg_b(X_i); | |
193 end; | |
194 x(X_i,end) = v; | |
195 end; | |
196 | |
197 Phi_x(:,end) = Phi(:,X_i)*x(X_i,end); | |
198 res = y - Phi_x(:,end); | |
199 if (~function_tau) % If tau is not a function handle... | |
200 if (optimizeTau) % Compute optimized tau | |
201 if (i > params.memory) | |
202 | |
203 % tau = argmin ||u - Phi*y_{i+1}|| | |
204 % = <res, Phi*(x_cur - x_prev)>/||Phi*(x_cur - x_prev)||^2 | |
205 | |
206 for k = 1:params.memory | |
207 Phi_diff = Phi_x(:,end) - Phi_x(:,k); | |
208 tau_candidates(k) = res'*Phi_diff/(Phi_diff'*Phi_diff); | |
209 y_candidates(:,k) = x(:,end) + tau_candidates(k)*(x(:,end) - x(:,k)); | |
210 residual_energy(k) = norm(y - Phi*y_candidates(:,k)); | |
211 end; | |
212 | |
213 [tmpArg,min_ind] = min(residual_energy); | |
214 y_cur = y_candidates(:,min_ind); | |
215 Y_i = find(ne(y_cur,0)); | |
216 else | |
217 Phi_diff = Phi_x(:,end) - Phi_x(:,end-1); | |
218 params.tau = res'*Phi_diff/(Phi_diff'*Phi_diff); | |
219 y_cur = x(:,end) + params.tau*(x(:,end) - x(:,end-1)); | |
220 Y_i = find(ne(y_cur,0)); | |
221 end; | |
222 else | |
223 if (i > params.memory) | |
224 for k = 1:params.memory | |
225 y_candidates(:,k) = x(:,end) + params.tau*(x(:,end) - x(:,k)); | |
226 residual_energy(k) = norm(y - Phi*y_candidates(:,k)); | |
227 end; | |
228 | |
229 [tmpArg,min_ind] = min(residual_energy); | |
230 y_cur = y_candidates(:,min_ind); | |
231 Y_i = find(ne(y_cur,0)); | |
232 else | |
233 y_cur = x(:,end) + params.tau*(x(:,end) - x(:,end-1)); | |
234 Y_i = find(ne(y_cur,0)); | |
235 end; | |
236 end; | |
237 else | |
238 if (i > params.memory) | |
239 for k = 1:params.memory | |
240 y_candidates(:,k) = x(:,end) + params.tau(i)*(x(:,end) - x(:,k)); | |
241 residual_energy(k) = norm(y - Phi*y_candidates(:,k))^2; | |
242 end; | |
243 | |
244 [tmpArg,min_ind] = min(residual_energy); | |
245 y_cur = y_candidates(:,min_ind); | |
246 Y_i = find(ne(y_cur,0)); | |
247 else | |
248 y_cur = x(:,end) + params.tau(i)*(x(:,end) - x(:,end-1)); | |
249 Y_i = find(ne(y_cur,0)); | |
250 end; | |
251 end; | |
252 | |
253 % Test stopping criterion | |
254 if (i > 1) && (norm(x(:,end) - x(:,end-1)) < params.tol*norm(x(:,end))) | |
255 break; | |
256 end; | |
257 i = i + 1; | |
258 end; | |
259 | |
260 x_hat = x(:,end); | |
261 numiter = i; | |
262 | |
263 if (i > params.ALPSiters) | |
264 x_path = x_path(:,1:numiter-1); | |
265 else | |
266 x_path = x_path(:,1:numiter); | |
267 end; |