Mercurial > hg > octave-nkf
annotate scripts/sparse/pcr.m @ 15536:2e8eb9ac43a5 stable rc-3-6-4-0
3.6.4-rc0 release candidate
* configure.ac (AC_INIT): Version is now 3.6.2-rc0.
(OCTAVE_RELEASE_DATE): Now 2012-05-11.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 17 Oct 2012 10:05:44 -0400 |
parents | 4d917a6a858b |
children | ce2b59a6d0e5 |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
14093
diff
changeset
|
1 ## Copyright (C) 2004-2012 Piotr Krzyzanowski |
5837 | 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 | |
7016 | 7 ## the Free Software Foundation; either version 3 of the License, or (at |
8 ## your option) any later version. | |
5837 | 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 | |
7016 | 16 ## along with Octave; see the file COPYING. If not, see |
17 ## <http://www.gnu.org/licenses/>. | |
5837 | 18 |
19 ## -*- texinfo -*- | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
20 ## @deftypefn {Function File} {@var{x} =} pcr (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{m}, @var{x0}, @dots{}) |
5837 | 21 ## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} pcr (@dots{}) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
22 ## |
14093
050bc580cb60
doc: Various docstring improvements before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
23 ## Solve the linear system of equations @code{@var{A} * @var{x} = @var{b}} |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
24 ## by means of the Preconditioned Conjugate Residuals iterative |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
25 ## method. The input arguments are |
5837 | 26 ## |
27 ## @itemize | |
28 ## @item | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
29 ## @var{A} can be either a square (preferably sparse) matrix or a |
5837 | 30 ## function handle, inline function or string containing the name |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
31 ## of a function which computes @code{@var{A} * @var{x}}. In principle |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
32 ## @var{A} should be symmetric and non-singular; if @code{pcr} |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
33 ## finds @var{A} to be numerically singular, you will get a warning |
5837 | 34 ## message and the @var{flag} output parameter will be set. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
35 ## |
5837 | 36 ## @item |
37 ## @var{b} is the right hand side vector. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
38 ## |
5837 | 39 ## @item |
40 ## @var{tol} is the required relative tolerance for the residual error, | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
41 ## @code{@var{b} - @var{A} * @var{x}}. The iteration stops if |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
42 ## @code{norm (@var{b} - @var{A} * @var{x}) <= |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
43 ## @var{tol} * norm (@var{b} - @var{A} * @var{x0})}. |
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
44 ## If @var{tol} is empty or is omitted, the function sets |
5837 | 45 ## @code{@var{tol} = 1e-6} by default. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
46 ## |
5837 | 47 ## @item |
48 ## @var{maxit} is the maximum allowable number of iterations; if | |
49 ## @code{[]} is supplied for @code{maxit}, or @code{pcr} has less | |
50 ## arguments, a default value equal to 20 is used. | |
51 ## | |
52 ## @item | |
5838 | 53 ## @var{m} is the (left) preconditioning matrix, so that the iteration is |
5837 | 54 ## (theoretically) equivalent to solving by @code{pcr} @code{@var{P} * |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
55 ## @var{x} = @var{m} \ @var{b}}, with @code{@var{P} = @var{m} \ @var{A}}. |
5837 | 56 ## Note that a proper choice of the preconditioner may dramatically |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
57 ## improve the overall performance of the method. Instead of matrix |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
58 ## @var{m}, the user may pass a function which returns the results of |
5838 | 59 ## applying the inverse of @var{m} to a vector (usually this is the |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
60 ## preferred way of using the preconditioner). If @code{[]} is supplied |
5838 | 61 ## for @var{m}, or @var{m} is omitted, no preconditioning is applied. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
62 ## |
5837 | 63 ## @item |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
64 ## @var{x0} is the initial guess. If @var{x0} is empty or omitted, the |
5837 | 65 ## function sets @var{x0} to a zero vector by default. |
66 ## @end itemize | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
67 ## |
5837 | 68 ## The arguments which follow @var{x0} are treated as parameters, and |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
69 ## passed in a proper way to any of the functions (@var{A} or @var{m}) |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
70 ## which are passed to @code{pcr}. See the examples below for further |
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
71 ## details. The output arguments are |
5837 | 72 ## |
73 ## @itemize | |
74 ## @item | |
75 ## @var{x} is the computed approximation to the solution of | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
76 ## @code{@var{A} * @var{x} = @var{b}}. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
77 ## |
5837 | 78 ## @item |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
79 ## @var{flag} reports on the convergence. @code{@var{flag} = 0} means |
5837 | 80 ## the solution converged and the tolerance criterion given by @var{tol} |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
81 ## is satisfied. @code{@var{flag} = 1} means that the @var{maxit} limit |
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
82 ## for the iteration count was reached. @code{@var{flag} = 3} reports t |
5837 | 83 ## @code{pcr} breakdown, see [1] for details. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
84 ## |
5837 | 85 ## @item |
86 ## @var{relres} is the ratio of the final residual to its initial value, | |
87 ## measured in the Euclidean norm. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
88 ## |
5837 | 89 ## @item |
90 ## @var{iter} is the actual number of iterations performed. | |
91 ## | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
92 ## @item |
5837 | 93 ## @var{resvec} describes the convergence history of the method, |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
94 ## so that @code{@var{resvec} (i)} contains the Euclidean norms of the |
7001 | 95 ## residual after the (@var{i}-1)-th iteration, @code{@var{i} = |
5838 | 96 ## 1,2, @dots{}, @var{iter}+1}. |
5837 | 97 ## @end itemize |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
98 ## |
5837 | 99 ## Let us consider a trivial problem with a diagonal matrix (we exploit the |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
100 ## sparsity of A) |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
101 ## |
5837 | 102 ## @example |
103 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
104 ## n = 10; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
105 ## A = sparse (diag (1:n)); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
106 ## b = rand (N, 1); |
5837 | 107 ## @end group |
108 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
109 ## |
5837 | 110 ## @sc{Example 1:} Simplest use of @code{pcr} |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
111 ## |
5837 | 112 ## @example |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
113 ## x = pcr (A, b) |
5837 | 114 ## @end example |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
115 ## |
5837 | 116 ## @sc{Example 2:} @code{pcr} with a function which computes |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
117 ## @code{@var{A} * @var{x}}. |
5837 | 118 ## |
119 ## @example | |
120 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
121 ## function y = apply_a (x) |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
122 ## y = [1:10]' .* x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
123 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
124 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
125 ## x = pcr ("apply_a", b) |
5837 | 126 ## @end group |
127 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
128 ## |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
129 ## @sc{Example 3:} Preconditioned iteration, with full diagnostics. The |
5837 | 130 ## preconditioner (quite strange, because even the original matrix |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
131 ## @var{A} is trivial) is defined as a function |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
132 ## |
5837 | 133 ## @example |
134 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
135 ## function y = apply_m (x) |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
136 ## k = floor (length (x) - 2); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
137 ## y = x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
138 ## y(1:k) = x(1:k) ./ [1:k]'; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
139 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
140 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
141 ## [x, flag, relres, iter, resvec] = ... |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
142 ## pcr (A, b, [], [], "apply_m") |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
143 ## semilogy ([1:iter+1], resvec); |
5837 | 144 ## @end group |
145 ## @end example | |
146 ## | |
147 ## @sc{Example 4:} Finally, a preconditioner which depends on a | |
148 ## parameter @var{k}. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
149 ## |
5837 | 150 ## @example |
151 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
152 ## function y = apply_m (x, varargin) |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
153 ## k = varargin@{1@}; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
154 ## y = x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
155 ## y(1:k) = x(1:k) ./ [1:k]'; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
156 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
157 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
158 ## [x, flag, relres, iter, resvec] = ... |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
159 ## pcr (A, b, [], [], "apply_m"', [], 3) |
5837 | 160 ## @end group |
161 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
162 ## |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
163 ## References: |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
164 ## |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
165 ## [1] W. Hackbusch, @cite{Iterative Solution of Large Sparse Systems of |
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
166 ## Equations}, section 9.5.4; Springer, 1994 |
5837 | 167 ## |
168 ## @seealso{sparse, pcg} | |
169 ## @end deftypefn | |
170 | |
5838 | 171 ## Author: Piotr Krzyzanowski <piotr.krzyzanowski@mimuw.edu.pl> |
5837 | 172 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
173 function [x, flag, relres, iter, resvec] = pcr (A, b, tol, maxit, m, x0, varargin) |
5837 | 174 |
175 breakdown = false; | |
176 | |
5838 | 177 if (nargin < 6 || isempty (x0)) |
178 x = zeros (size (b)); | |
5837 | 179 else |
180 x = x0; | |
181 endif | |
182 | |
183 if (nargin < 5) | |
8507 | 184 m = []; |
5837 | 185 endif |
186 | |
5838 | 187 if (nargin < 4 || isempty (maxit)) |
5837 | 188 maxit = 20; |
189 endif | |
190 | |
5838 | 191 maxit += 2; |
5837 | 192 |
5838 | 193 if (nargin < 3 || isempty (tol)) |
5837 | 194 tol = 1e-6; |
195 endif | |
196 | |
197 if (nargin < 2) | |
5838 | 198 print_usage (); |
5837 | 199 endif |
200 | |
201 ## init | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
202 if (isnumeric (A)) # is A a matrix? |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
203 r = b - A*x; |
10549 | 204 else # then A should be a function! |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
205 r = b - feval (A, x, varargin{:}); |
5837 | 206 endif |
207 | |
10549 | 208 if (isnumeric (m)) # is M a matrix? |
209 if (isempty (m)) # if M is empty, use no precond | |
5837 | 210 p = r; |
10549 | 211 else # otherwise, apply the precond |
8507 | 212 p = m \ r; |
5837 | 213 endif |
10549 | 214 else # then M should be a function! |
8507 | 215 p = feval (m, r, varargin{:}); |
5837 | 216 endif |
217 | |
218 iter = 2; | |
219 | |
220 b_bot_old = 1; | |
5838 | 221 q_old = p_old = s_old = zeros (size (x)); |
5837 | 222 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
223 if (isnumeric (A)) # is A a matrix? |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
224 q = A * p; |
10549 | 225 else # then A should be a function! |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
226 q = feval (A, p, varargin{:}); |
5837 | 227 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
228 |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
229 resvec(1) = abs (norm (r)); |
5837 | 230 |
231 ## iteration | |
5838 | 232 while (resvec(iter-1) > tol*resvec(1) && iter < maxit) |
233 | |
10549 | 234 if (isnumeric (m)) # is M a matrix? |
235 if (isempty (m)) # if M is empty, use no precond | |
236 s = q; | |
237 else # otherwise, apply the precond | |
238 s = m \ q; | |
5837 | 239 endif |
10549 | 240 else # then M should be a function! |
8507 | 241 s = feval (m, q, varargin{:}); |
5837 | 242 endif |
5838 | 243 b_top = r' * s; |
244 b_bot = q' * s; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
245 |
5837 | 246 if (b_bot == 0.0) |
247 breakdown = true; | |
248 break; | |
249 endif | |
5838 | 250 lambda = b_top / b_bot; |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
251 |
5838 | 252 x += lambda*p; |
253 r -= lambda*q; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
254 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
255 if (isnumeric(A)) # is A a matrix? |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
256 t = A*s; |
10549 | 257 else # then A should be a function! |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
258 t = feval (A, s, varargin{:}); |
5837 | 259 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
260 |
5838 | 261 alpha0 = (t'*s) / b_bot; |
262 alpha1 = (t'*s_old) / b_bot_old; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
263 |
5838 | 264 p_temp = p; |
265 q_temp = q; | |
266 | |
5837 | 267 p = s - alpha0*p - alpha1*p_old; |
268 q = t - alpha0*q - alpha1*q_old; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
269 |
5837 | 270 s_old = s; |
271 p_old = p_temp; | |
272 q_old = q_temp; | |
273 b_bot_old = b_bot; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
274 |
5838 | 275 resvec(iter) = abs (norm (r)); |
276 iter++; | |
5837 | 277 endwhile |
278 | |
279 flag = 0; | |
5838 | 280 relres = resvec(iter-1) ./ resvec(1); |
281 iter -= 2; | |
282 if (iter >= maxit-2) | |
5837 | 283 flag = 1; |
284 if (nargout < 2) | |
6498 | 285 warning ("pcr: maximum number of iterations (%d) reached\n", iter); |
286 warning ("the initial residual norm was reduced %g times.\n", 1.0/relres); | |
5837 | 287 endif |
5838 | 288 elseif (nargout < 2 && ! breakdown) |
6498 | 289 fprintf (stderr, "pcr: converged in %d iterations. \n", iter); |
290 fprintf (stderr, "the initial residual norm was reduced %g times.\n", | |
10549 | 291 1.0 / relres); |
5837 | 292 endif |
5838 | 293 |
5837 | 294 if (breakdown) |
295 flag = 3; | |
296 if (nargout < 2) | |
7001 | 297 warning ("pcr: breakdown occurred:\n"); |
6498 | 298 warning ("system matrix singular or preconditioner indefinite?\n"); |
5837 | 299 endif |
300 endif | |
301 | |
302 endfunction | |
303 | |
304 %!demo | |
305 %! | |
10549 | 306 %! # Simplest usage of PCR (see also 'help pcr') |
5837 | 307 %! |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
308 %! N = 20; |
10549 | 309 %! A = diag(linspace(-3.1,3,N)); b = rand(N,1); y = A\b; #y is the true solution |
310 %! x = pcr(A,b); | |
311 %! printf('The solution relative error is %g\n', norm(x-y)/norm(y)); | |
5837 | 312 %! |
10549 | 313 %! # You shouldn't be afraid if PCR issues some warning messages in this |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
314 %! # example: watch out in the second example, why it takes N iterations |
10549 | 315 %! # of PCR to converge to (a very accurate, by the way) solution |
5837 | 316 %!demo |
317 %! | |
10549 | 318 %! # Full output from PCR |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
319 %! # We use this output to plot the convergence history |
5837 | 320 %! |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
321 %! N = 20; |
10549 | 322 %! A = diag(linspace(-3.1,30,N)); b = rand(N,1); X = A\b; #X is the true solution |
323 %! [x, flag, relres, iter, resvec] = pcr(A,b); | |
324 %! printf('The solution relative error is %g\n', norm(x-X)/norm(X)); | |
325 %! title('Convergence history'); xlabel('Iteration'); ylabel('log(||b-Ax||/||b||)'); | |
326 %! semilogy([0:iter],resvec/resvec(1),'o-g;relative residual;'); | |
5837 | 327 %!demo |
328 %! | |
10549 | 329 %! # Full output from PCR |
330 %! # We use indefinite matrix based on the Hilbert matrix, with one | |
331 %! # strongly negative eigenvalue | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
332 %! # Hilbert matrix is extremely ill conditioned, so is ours, |
10549 | 333 %! # and that's why PCR WILL have problems |
5837 | 334 %! |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
335 %! N = 10; |
10549 | 336 %! A = hilb(N); A(1,1)=-A(1,1); b = rand(N,1); X = A\b; #X is the true solution |
337 %! printf('Condition number of A is %g\n', cond(A)); | |
338 %! [x, flag, relres, iter, resvec] = pcr(A,b,[],200); | |
339 %! if (flag == 3) | |
340 %! printf('PCR breakdown. System matrix is [close to] singular\n'); | |
341 %! end | |
342 %! title('Convergence history'); xlabel('Iteration'); ylabel('log(||b-Ax||)'); | |
343 %! semilogy([0:iter],resvec,'o-g;absolute residual;'); | |
5837 | 344 %!demo |
345 %! | |
10549 | 346 %! # Full output from PCR |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
347 %! # We use an indefinite matrix based on the 1-D Laplacian matrix for A, |
10549 | 348 %! # and here we have cond(A) = O(N^2) |
349 %! # That's the reason we need some preconditioner; here we take | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
350 %! # a very simple and not powerful Jacobi preconditioner, |
10549 | 351 %! # which is the diagonal of A |
5837 | 352 %! |
10549 | 353 %! # Note that we use here indefinite preconditioners! |
5837 | 354 %! |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
355 %! N = 100; |
10549 | 356 %! A = zeros(N,N); |
357 %! for i=1:N-1 # form 1-D Laplacian matrix | |
358 %! A(i:i+1,i:i+1) = [2 -1; -1 2]; | |
359 %! endfor | |
360 %! A = [A, zeros(size(A)); zeros(size(A)), -A]; | |
361 %! b = rand(2*N,1); X = A\b; #X is the true solution | |
362 %! maxit = 80; | |
363 %! printf('System condition number is %g\n',cond(A)); | |
364 %! # No preconditioner: the convergence is very slow! | |
5837 | 365 %! |
10549 | 366 %! [x, flag, relres, iter, resvec] = pcr(A,b,[],maxit); |
367 %! title('Convergence history'); xlabel('Iteration'); ylabel('log(||b-Ax||)'); | |
368 %! semilogy([0:iter],resvec,'o-g;NO preconditioning: absolute residual;'); | |
5837 | 369 %! |
10549 | 370 %! pause(1); |
371 %! # Test Jacobi preconditioner: it will not help much!!! | |
5837 | 372 %! |
10549 | 373 %! M = diag(diag(A)); # Jacobi preconditioner |
374 %! [x, flag, relres, iter, resvec] = pcr(A,b,[],maxit,M); | |
375 %! hold on; | |
376 %! semilogy([0:iter],resvec,'o-r;JACOBI preconditioner: absolute residual;'); | |
5837 | 377 %! |
10549 | 378 %! pause(1); |
379 %! # Test nonoverlapping block Jacobi preconditioner: this one should give | |
380 %! # some convergence speedup! | |
5837 | 381 %! |
10549 | 382 %! M = zeros(N,N);k=4; |
383 %! for i=1:k:N # get k x k diagonal blocks of A | |
384 %! M(i:i+k-1,i:i+k-1) = A(i:i+k-1,i:i+k-1); | |
385 %! endfor | |
386 %! M = [M, zeros(size(M)); zeros(size(M)), -M]; | |
387 %! [x, flag, relres, iter, resvec] = pcr(A,b,[],maxit,M); | |
388 %! semilogy([0:iter],resvec,'o-b;BLOCK JACOBI preconditioner: absolute residual;'); | |
389 %! hold off; | |
5837 | 390 %!test |
391 %! | |
10549 | 392 %! #solve small indefinite diagonal system |
5837 | 393 %! |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
394 %! N = 10; |
10549 | 395 %! A = diag(linspace(-10.1,10,N)); b = ones(N,1); X = A\b; #X is the true solution |
396 %! [x, flag] = pcr(A,b,[],N+1); | |
397 %! assert(norm(x-X)/norm(X)<1e-10); | |
398 %! assert(flag,0); | |
5837 | 399 %! |
400 %!test | |
401 %! | |
10549 | 402 %! #solve tridiagonal system, do not converge in default 20 iterations |
403 %! #should perform max allowable default number of iterations | |
5837 | 404 %! |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
405 %! N = 100; |
10549 | 406 %! A = zeros(N,N); |
407 %! for i=1:N-1 # form 1-D Laplacian matrix | |
408 %! A(i:i+1,i:i+1) = [2 -1; -1 2]; | |
409 %! endfor | |
410 %! b = ones(N,1); X = A\b; #X is the true solution | |
411 %! [x, flag, relres, iter, resvec] = pcr(A,b,1e-12); | |
412 %! assert(flag,1); | |
413 %! assert(relres>0.6); | |
414 %! assert(iter,20); | |
5837 | 415 %! |
416 %!test | |
417 %! | |
10549 | 418 %! #solve tridiagonal system with 'prefect' preconditioner |
419 %! #converges in one iteration | |
5837 | 420 %! |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
421 %! N = 100; |
10549 | 422 %! A = zeros(N,N); |
423 %! for i=1:N-1 # form 1-D Laplacian matrix | |
424 %! A(i:i+1,i:i+1) = [2 -1; -1 2]; | |
425 %! endfor | |
426 %! b = ones(N,1); X = A\b; #X is the true solution | |
427 %! [x, flag, relres, iter] = pcr(A,b,[],[],A,b); | |
428 %! assert(norm(x-X)/norm(X)<1e-6); | |
429 %! assert(relres<1e-6); | |
430 %! assert(flag,0); | |
431 %! assert(iter,1); #should converge in one iteration | |
5837 | 432 %! |