5837
|
1 ## Copyright (C) 2004 Piotr Krzyzanowski <piotr.krzyzanowski@mimuw.edu.pl> |
|
2 ## |
|
3 ## This file is part of Octave. |
|
4 ## |
|
5 ## Octave is free software; you can redistribute it and/or modify it |
|
6 ## under the terms of the GNU General Public License as published by |
|
7 ## the Free Software Foundation; either version 2, or (at your option) |
|
8 ## any later version. |
|
9 ## |
|
10 ## Octave is distributed in the hope that it will be useful, but |
|
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
|
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
|
13 ## General Public License for more details. |
|
14 ## |
|
15 ## You should have received a copy of the GNU General Public License |
|
16 ## along with Octave; see the file COPYING. If not, write to the Free |
|
17 ## Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA |
|
18 ## 02110-1301, USA. |
|
19 |
|
20 ## -*- texinfo -*- |
|
21 ## @deftypefn {Function File} {@var{x} =} pcg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M}, @var{x0}, @dots{}) |
|
22 ## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}, @var{eigest}] =} pcg (@dots{}) |
|
23 ## |
|
24 ## Solves the linear system of equations @code{@var{A} * @var{x} = |
|
25 ## @var{b}} by means of the Preconditioned Conjugate Gradient iterative |
|
26 ## method. The input arguments are |
|
27 ## |
|
28 ## @itemize |
|
29 ## @item |
|
30 ## @var{A} can be either a square (preferably sparse) matrix or a |
|
31 ## function handle, inline function or string containing the name |
|
32 ## of a function which computes @code{@var{A} * @var{x}}. In principle |
|
33 ## @var{A} should be symmetric and positive definite; if @code{pcg} |
|
34 ## finds @var{A} to not be positive definite, you will get a warning |
|
35 ## message and the @var{flag} output parameter will be set. |
|
36 ## |
|
37 ## @item |
|
38 ## @var{b} is the right hand side vector. |
|
39 ## |
|
40 ## @item |
|
41 ## @var{tol} is the required relative tolerance for the residual error, |
|
42 ## @code{@var{b} - @var{A} * @var{x}}. The iteration stops if @code{norm |
|
43 ## (@var{b} - @var{A} * @var{x}) <= @var{tol} * norm (@var{b} - @var{A} * |
|
44 ## @var{x0})}. If @var{tol} is empty or is omitted, the function sets |
|
45 ## @code{@var{tol} = 1e-6} by default. |
|
46 ## |
|
47 ## @item |
|
48 ## @var{maxit} is the maximum allowable number of iterations; if |
|
49 ## @code{[]} is supplied for @code{maxit}, or @code{pcg} has less |
|
50 ## arguments, a default value equal to 20 is used. |
|
51 ## |
|
52 ## @item |
|
53 ## @var{M} is the (left) preconditioning matrix, so that the iteration is |
|
54 ## (theoretically) equivalent to solving by @code{pcg} @code{@var{P} * |
|
55 ## @var{x} = @var{M} \ @var{b}}, with @code{@var{P} = @var{M} \ @var{A}}. |
|
56 ## Note that a proper choice of the preconditioner may dramatically |
|
57 ## improve the overall performance of the method. Instead of matrix |
|
58 ## @var{M}, the user may pass a function which returns the results of |
|
59 ## applying the inverse of @var{M} to a vector (usually this is the |
|
60 ## preferred way of using the preconditioner). If @code{[]} is supplied |
|
61 ## for @var{M}, or @var{M} is omitted, no preconditioning is applied. |
|
62 ## |
|
63 ## @item |
|
64 ## @var{x0} is the initial guess. If @var{x0} is empty or omitted, the |
|
65 ## function sets @var{x0} to a zero vector by default. |
|
66 ## @end itemize |
|
67 ## |
|
68 ## The arguments which follow @var{x0} are treated as parameters, and |
|
69 ## passed in a proper way to any of the functions (@var{A} or @var{M}) |
|
70 ## which are passed to @code{pcg}. See the examples below for further |
|
71 ## details. The output arguments are |
|
72 ## |
|
73 ## @itemize |
|
74 ## @item |
|
75 ## @var{x} is the computed approximation to the solution of |
|
76 ## @code{@var{A} * @var{x} = @var{b}}. |
|
77 ## |
|
78 ## @item |
|
79 ## @var{flag} reports on the convergence. @code{@var{flag} = 0} means |
|
80 ## the solution converged and the tolerance criterion given by @var{tol} |
|
81 ## is satisfied. @code{@var{flag} = 1} means that the @var{maxit} limit |
|
82 ## for the iteration count was reached. @code{@var{flag} = 3} reports that |
|
83 ## the (preconditioned) matrix was found not positive definite. |
|
84 ## |
|
85 ## @item |
|
86 ## @var{relres} is the ratio of the final residual to its initial value, |
|
87 ## measured in the Euclidean norm. |
|
88 ## |
|
89 ## @item |
|
90 ## @var{iter} is the actual number of iterations performed. |
|
91 ## |
|
92 ## @item |
|
93 ## @var{resvec} describes the convergence history of the method. |
|
94 ## @code{@var{resvec} (i,1)} is the Euclidean norm of the residual, and |
|
95 ## @code{@var{resvec} (i,2)} is the preconditioned residual norm, |
|
96 ## after the (@var{i}-1)-th iteration, @code{@var{i} = |
|
97 ## 1,2,...@var{iter}+1}. The preconditioned residual norm is defined as |
|
98 ## @code{norm (@var{r}) ^ 2 = @var{r}' * (@var{M} \ @var{r})} where |
|
99 ## @code{@var{r} = @var{b} - @var{A} * @var{x}}, see also the |
|
100 ## description of @var{M}. If @var{eigest} is not required, only |
|
101 ## @code{@var{resvec} (:,1)} is returned. |
|
102 ## |
|
103 ## @item |
|
104 ## @var{eigest} returns the estimate for the smallest @code{@var{eigest} |
|
105 ## (1)} and largest @code{@var{eigest} (2)} eigenvalues of the |
|
106 ## preconditioned matrix @code{@var{P} = @var{M} \ @var{A}}. In |
|
107 ## particular, if no preconditioning is used, the extimates for the |
|
108 ## extreme eigenvalues of @var{A} are returned. @code{@var{eigest} (1)} |
|
109 ## is an overestimate and @code{@var{eigest} (2)} is an underestimate, |
|
110 ## so that @code{@var{eigest} (2) / @var{eigest} (1)} is a lower bound |
|
111 ## for @code{cond (@var{P}, 2)}, which nevertheless in the limit should |
|
112 ## theoretically be equal to the actual value of the condition number. |
|
113 ## The method which computes @var{eigest} works only for symmetric positive |
|
114 ## definite @var{A} and @var{M}, and the user is responsible for |
|
115 ## verifying this assumption. |
|
116 ## @end itemize |
|
117 ## |
|
118 ## Let us consider a trivial problem with a diagonal matrix (we exploit the |
|
119 ## sparsity of A) |
|
120 ## |
|
121 ## @example |
|
122 ## @group |
|
123 ## N = 10; |
|
124 ## A = diag([1:N]); A = sparse(A); |
|
125 ## b = rand(N,1); |
|
126 ## @end group |
|
127 ## @end example |
|
128 ## |
|
129 ## @sc{Example 1:} Simplest use of @code{pcg} |
|
130 ## |
|
131 ## @example |
|
132 ## x = pcg(A,b) |
|
133 ## @end example |
|
134 ## |
|
135 ## @sc{Example 2:} @code{pcg} with a function which computes |
|
136 ## @code{@var{A} * @var{x}} |
|
137 ## |
|
138 ## @example |
|
139 ## @group |
|
140 ## function y = applyA(x) |
|
141 ## y = [1:N]'.*x; |
|
142 ## endfunction |
|
143 ## |
|
144 ## x = pcg('applyA',b) |
|
145 ## @end group |
|
146 ## @end example |
|
147 ## |
|
148 ## @sc{Example 3:} Preconditioned iteration, with full diagnostics. The |
|
149 ## preconditioner (quite strange, because even the original matrix |
|
150 ## @var{A} is trivial) is defined as a function |
|
151 ## |
|
152 ## @example |
|
153 ## @group |
|
154 ## function y = applyM(x) |
|
155 ## K = floor(length(x)-2); |
|
156 ## y = x; |
|
157 ## y(1:K) = x(1:K)./[1:K]'; |
|
158 ## endfunction |
|
159 ## |
|
160 ## [x, flag, relres, iter, resvec, eigest] = pcg(A,b,[],[],'applyM') |
|
161 ## semilogy([1:iter+1], resvec); |
|
162 ## @end group |
|
163 ## @end example |
|
164 ## |
|
165 ## @sc{Example 4:} Finally, a preconditioner which depends on a |
|
166 ## parameter @var{k}. |
|
167 ## |
|
168 ## @example |
|
169 ## @group |
|
170 ## function y = applyM(x, varargin) |
|
171 ## K = varargin@{1@}; |
|
172 ## y = x; y(1:K) = x(1:K)./[1:K]'; |
|
173 ## endfuntion |
|
174 ## |
|
175 ## [x, flag, relres, iter, resvec, eigest] = ... |
|
176 ## pcg(A,b,[],[],'applyM',[],3) |
|
177 ## @end group |
|
178 ## @end example |
|
179 ## |
|
180 ## @sc{References} |
|
181 ## |
|
182 ## [1] C.T.Kelley, 'Iterative methods for linear and nonlinear equations', |
|
183 ## SIAM, 1995 (the base PCG algorithm) |
|
184 ## |
|
185 ## [2] Y.Saad, 'Iterative methods for sparse linear systems', PWS 1996 |
|
186 ## (condition number estimate from PCG) Revised version of this book is |
|
187 ## available online at http://www-users.cs.umn.edu/~saad/books.html |
|
188 ## |
|
189 ## |
|
190 ## @seealso{sparse, pcr} |
|
191 ## @end deftypefn |
|
192 |
|
193 ## REVISION HISTORY |
|
194 ## |
|
195 ## 2004-05-21, Piotr Krzyzanowski: |
|
196 ## Added 4 demos and 4 tests |
|
197 ## |
|
198 ## 2004-05-18, Piotr Krzyzanowski: |
|
199 ## Warnings use warning() function now |
|
200 ## |
|
201 ## 2004-04-29, Piotr Krzyzanowski: |
|
202 ## Added more warning messages when FLAG is not required |
|
203 ## |
|
204 ## 2004-04-28, Piotr Krzyzanowski: |
|
205 ## When eigest is required, resvec returns both the Euclidean and the |
|
206 ## preconditioned residual norm convergence history |
|
207 ## |
|
208 ## 2004-04-20, Piotr Krzyzanowski: |
|
209 ## Corrected eigenvalue estimator. Changed the tridiagonal matrix |
|
210 ## eigenvalue solver to regular eig |
|
211 ## |
|
212 |
|
213 function [x, flag, relres, iter, resvec, eigest] = ... |
|
214 pcg( A, b, tol, maxit, M, x0, varargin ) |
|
215 |
|
216 if ((nargin < 6) || isempty(x0)) |
|
217 x = zeros(size(b)); |
|
218 else |
|
219 x = x0; |
|
220 endif |
|
221 |
|
222 if (nargin < 5) |
|
223 M = []; |
|
224 endif |
|
225 |
|
226 if ((nargin < 4) || isempty(maxit)) |
|
227 maxit = min(size(b,1),20); |
|
228 endif |
|
229 |
|
230 maxit = maxit + 2; |
|
231 |
|
232 if ((nargin < 3) || isempty(tol)) |
|
233 tol = 1e-6; |
|
234 endif |
|
235 |
|
236 preconditioned_residual_out = false; |
|
237 if (nargout > 5) |
|
238 T = zeros(maxit,maxit); |
|
239 preconditioned_residual_out = true; |
|
240 endif |
|
241 |
|
242 matrix_positive_definite = true; # assume A is positive definite |
|
243 |
|
244 p = zeros(size(b)); |
|
245 oldtau = 1; |
|
246 if (isnumeric(A)) # is A a matrix? |
|
247 r = b - A*x; |
|
248 else # then A should be a function! |
|
249 r = b - feval(A,x,varargin{:}); |
|
250 endif |
|
251 |
|
252 resvec(1,1) = norm(r); |
|
253 alpha = 1; |
|
254 iter = 2; |
|
255 |
|
256 while ((resvec(iter-1,1) > tol*resvec(1,1)) && (iter < maxit)) |
|
257 if (isnumeric(M)) # is M a matrix? |
|
258 if isempty(M) # if M is empty, use no precond |
|
259 z = r; |
|
260 else # otherwise, apply the precond |
|
261 z = M \ r; |
|
262 endif |
|
263 else # then M should be a function! |
|
264 z = feval(M,r,varargin{:}); |
|
265 endif |
|
266 tau = z'*r; |
|
267 resvec(iter-1,2) = sqrt(tau); |
|
268 beta = tau/oldtau; |
|
269 oldtau = tau; |
|
270 p = z + beta*p; |
|
271 if (isnumeric(A)) # is A a matrix? |
|
272 w = A*p; |
|
273 else # then A should be a function! |
|
274 w = feval(A,p,varargin{:}); |
|
275 endif |
|
276 oldalpha = alpha; # needed only for eigest |
|
277 alpha = tau/(p'*w); |
|
278 if (alpha <= 0.0) # negative matrix? |
|
279 matrix_positive_definite = false; |
|
280 endif |
|
281 x = x + alpha*p; |
|
282 r = r - alpha*w; |
|
283 if ((nargout > 5) && (iter > 2)) |
|
284 T(iter-1:iter, iter-1:iter) = T(iter-1:iter, iter-1:iter) + ... |
|
285 [1 sqrt(beta); sqrt(beta) beta]./oldalpha; |
|
286 ## EVS = eig(T(2:iter-1,2:iter-1)); |
|
287 ## fprintf(stderr,"PCG condest: %g (iteration: %d)\n", max(EVS)/min(EVS),iter); |
|
288 endif |
|
289 resvec(iter,1) = norm(r); |
|
290 iter = iter + 1; |
|
291 endwhile |
|
292 |
|
293 if (nargout > 5) |
|
294 if (matrix_positive_definite ) |
|
295 if (iter > 3) |
|
296 T = T(2:iter-2,2:iter-2); |
|
297 l = eig(T); |
|
298 eigest = [min(l) max(l)]; |
|
299 ## fprintf(stderr, "PCG condest: %g\n",eigest(2)/eigest(1)); |
|
300 else |
|
301 eigest = [NaN NaN]; |
|
302 warning("PCG: eigenvalue estimate failed: iteration converged too fast."); |
|
303 endif |
|
304 else |
|
305 eigest = [NaN NaN]; |
|
306 endif |
|
307 |
|
308 ## apply the preconditioner once more and finish with the precond |
|
309 ## residual |
|
310 if (isnumeric(M)) # is M a matrix? |
|
311 if isempty(M) # if M is empty, use no precond |
|
312 z = r; |
|
313 else # otherwise, apply the precond |
|
314 z = M\r; |
|
315 endif |
|
316 else # then M should be a function! |
|
317 z = feval(M,r,varargin{:}); |
|
318 endif |
|
319 resvec(iter-1,2) = sqrt(r'*z); |
|
320 else |
|
321 resvec = resvec(:,1); |
|
322 endif |
|
323 |
|
324 flag = 0; |
|
325 relres = resvec(iter-1,1)./resvec(1,1); |
|
326 iter = iter - 2; |
|
327 if (iter >= (maxit-2)) |
|
328 flag = 1; |
|
329 if (nargout < 2) |
|
330 warning("PCG: maximum number of iterations (%d) reached\n", iter); |
|
331 warning("The initial residual norm was reduced %g times.\n", 1.0/relres); |
|
332 endif |
|
333 else |
|
334 if (nargout < 2) |
|
335 fprintf(stderr, "PCG: converged in %d iterations. ", iter); |
|
336 fprintf(stderr, "The initial residual norm was reduced %g times.\n",... |
|
337 1.0/relres); |
|
338 endif |
|
339 endif |
|
340 |
|
341 if (!matrix_positive_definite) |
|
342 flag = 3; |
|
343 if (nargout < 2) |
|
344 warning("PCG: matrix not positive definite?\n"); |
|
345 endif |
|
346 endif |
|
347 endfunction |
|
348 |
|
349 %!demo |
|
350 %! |
|
351 %! # Simplest usage of pcg (see also 'help pcg') |
|
352 %! |
|
353 %! N = 10; |
|
354 %! A = diag([1:N]); b = rand(N,1); y = A\b; #y is the true solution |
|
355 %! x = pcg(A,b); |
|
356 %! printf('The solution relative error is %g\n', norm(x-y)/norm(y)); |
|
357 %! |
|
358 %! # You shouldn't be afraid if pcg issues some warning messages in this |
|
359 %! # example: watch out in the second example, why it takes N iterations |
|
360 %! # of pcg to converge to (a very accurate, by the way) solution |
|
361 %!demo |
|
362 %! |
|
363 %! # Full output from pcg, except for the eigenvalue estimates |
|
364 %! # We use this output to plot the convergence history |
|
365 %! |
|
366 %! N = 10; |
|
367 %! A = diag([1:N]); b = rand(N,1); X = A\b; #X is the true solution |
|
368 %! [x, flag, relres, iter, resvec] = pcg(A,b); |
|
369 %! printf('The solution relative error is %g\n', norm(x-X)/norm(X)); |
|
370 %! title('Convergence history'); xlabel('Iteration'); ylabel('log(||b-Ax||/||b||)'); |
|
371 %! semilogy([0:iter],resvec/resvec(1),'o-g;relative residual;'); |
|
372 %!demo |
|
373 %! |
|
374 %! # Full output from pcg, including the eigenvalue estimates |
|
375 %! # Hilbert matrix is extremely ill conditioned, so pcg WILL have problems |
|
376 %! |
|
377 %! N = 10; |
|
378 %! A = hilb(N); b = rand(N,1); X = A\b; #X is the true solution |
|
379 %! [x, flag, relres, iter, resvec, eigest] = pcg(A,b,[],200); |
|
380 %! printf('The solution relative error is %g\n', norm(x-X)/norm(X)); |
|
381 %! printf('Condition number estimate is %g\n', eigest(2)/eigest(1)); |
|
382 %! printf('Actual condition number is %g\n', cond(A)); |
|
383 %! title('Convergence history'); xlabel('Iteration'); ylabel('log(||b-Ax||)'); |
|
384 %! semilogy([0:iter],resvec,['o-g;absolute residual;';'+-r;absolute preconditioned residual;']); |
|
385 %!demo |
|
386 %! |
|
387 %! # Full output from pcg, including the eigenvalue estimates |
|
388 %! # We use the 1-D Laplacian matrix for A, and cond(A) = O(N^2) |
|
389 %! # and that's the reasone we need some preconditioner; here we take |
|
390 %! # a very simple and not powerful Jacobi preconditioner, |
|
391 %! # which is the diagonal of A |
|
392 %! |
|
393 %! N = 100; |
|
394 %! A = zeros(N,N); |
|
395 %! for i=1:N-1 # form 1-D Laplacian matrix |
|
396 %! A(i:i+1,i:i+1) = [2 -1; -1 2]; |
|
397 %! endfor |
|
398 %! b = rand(N,1); X = A\b; #X is the true solution |
|
399 %! maxit = 80; |
|
400 %! printf('System condition number is %g\n',cond(A)); |
|
401 %! # No preconditioner: the convergence is very slow! |
|
402 %! |
|
403 %! [x, flag, relres, iter, resvec, eigest] = pcg(A,b,[],maxit); |
|
404 %! printf('System condition number estimate is %g\n',eigest(2)/eigest(1)); |
|
405 %! title('Convergence history'); xlabel('Iteration'); ylabel('log(||b-Ax||)'); |
|
406 %! semilogy([0:iter],resvec(:,1),'o-g;NO preconditioning: absolute residual;'); |
|
407 %! |
|
408 %! pause(1); |
|
409 %! # Test Jacobi preconditioner: it will not help much!!! |
|
410 %! |
|
411 %! M = diag(diag(A)); # Jacobi preconditioner |
|
412 %! [x, flag, relres, iter, resvec, eigest] = pcg(A,b,[],maxit,M); |
|
413 %! printf('JACOBI preconditioned system condition number estimate is %g\n',eigest(2)/eigest(1)); |
|
414 %! hold on; |
|
415 %! semilogy([0:iter],resvec(:,1),'o-r;JACOBI preconditioner: absolute residual;'); |
|
416 %! |
|
417 %! pause(1); |
|
418 %! # Test nonoverlapping block Jacobi preconditioner: it will help much! |
|
419 %! |
|
420 %! M = zeros(N,N);k=4 |
|
421 %! for i=1:k:N # form 1-D Laplacian matrix |
|
422 %! M(i:i+k-1,i:i+k-1) = A(i:i+k-1,i:i+k-1); |
|
423 %! endfor |
|
424 %! [x, flag, relres, iter, resvec, eigest] = pcg(A,b,[],maxit,M); |
|
425 %! printf('BLOCK JACOBI preconditioned system condition number estimate is %g\n',eigest(2)/eigest(1)); |
|
426 %! semilogy([0:iter],resvec(:,1),'o-b;BLOCK JACOBI preconditioner: absolute residual;'); |
|
427 %! hold off; |
|
428 %!test |
|
429 %! |
|
430 %! #solve small diagonal system |
|
431 %! |
|
432 %! N = 10; |
|
433 %! A = diag([1:N]); b = rand(N,1); X = A\b; #X is the true solution |
|
434 %! [x, flag] = pcg(A,b,[],N+1); |
|
435 %! assert(norm(x-X)/norm(X),0,1e-10); |
|
436 %! assert(flag,0); |
|
437 %! |
|
438 %!test |
|
439 %! |
|
440 %! #solve small indefinite diagonal system |
|
441 %! #despite A is indefinite, the iteration continues and converges |
|
442 %! #indefiniteness of A is detected |
|
443 %! |
|
444 %! N = 10; |
|
445 %! A = diag([1:N].*(-ones(1,N).^2)); b = rand(N,1); X = A\b; #X is the true solution |
|
446 %! [x, flag] = pcg(A,b,[],N+1); |
|
447 %! assert(norm(x-X)/norm(X),0,1e-10); |
|
448 %! assert(flag,3); |
|
449 %! |
|
450 %!test |
|
451 %! |
|
452 %! #solve tridiagonal system, do not converge in default 20 iterations |
|
453 %! |
|
454 %! N = 100; |
|
455 %! A = zeros(N,N); |
|
456 %! for i=1:N-1 # form 1-D Laplacian matrix |
|
457 %! A(i:i+1,i:i+1) = [2 -1; -1 2]; |
|
458 %! endfor |
|
459 %! b = ones(N,1); X = A\b; #X is the true solution |
|
460 %! [x, flag, relres, iter, resvec, eigest] = pcg(A,b,1e-12); |
|
461 %! assert(flag); |
|
462 %! assert(relres>1.0); |
|
463 %! assert(iter,20); #should perform max allowable default number of iterations |
|
464 %! |
|
465 %!test |
|
466 %! |
|
467 %! #solve tridiagonal system with 'prefect' preconditioner |
|
468 %! #converges in one iteration, so the eigest does not work |
|
469 %! #and issues a warning |
|
470 %! |
|
471 %! N = 100; |
|
472 %! A = zeros(N,N); |
|
473 %! for i=1:N-1 # form 1-D Laplacian matrix |
|
474 %! A(i:i+1,i:i+1) = [2 -1; -1 2]; |
|
475 %! endfor |
|
476 %! b = ones(N,1); X = A\b; #X is the true solution |
|
477 %! [x, flag, relres, iter, resvec, eigest] = pcg(A,b,[],[],A,b); |
|
478 %! assert(norm(x-X)/norm(X),0,1e-6); |
|
479 %! assert(flag,0); |
|
480 %! assert(iter,1); #should converge in one iteration |
|
481 %! assert(isnan(eigest),isnan([NaN NaN])); |
|
482 %! |