Mercurial > hg > octave-lyh
annotate scripts/sparse/pcg.m @ 17535:c12c688a35ed default tip lyh
Fix warnings
author | LYH <lyh.kernel@gmail.com> |
---|---|
date | Fri, 27 Sep 2013 17:43:27 +0800 |
parents | b81b9d079515 |
children |
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:
10846
diff
changeset
|
20 ## @deftypefn {Function File} {@var{x} =} pcg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{m1}, @var{m2}, @var{x0}, @dots{}) |
5837 | 21 ## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}, @var{eigest}] =} pcg (@dots{}) |
22 ## | |
14093
050bc580cb60
doc: Various docstring improvements before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
11588
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 Gradient 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:
10846
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:
10846
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:
10846
diff
changeset
|
32 ## @var{A} should be symmetric and positive definite; if @code{pcg} |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
33 ## finds @var{A} to not be positive definite, 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, | |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
41 ## @code{@var{b} - @var{A} * @var{x}}. The iteration stops if |
11587
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{pcg} has less | |
50 ## arguments, a default value equal to 20 is used. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
51 ## |
5837 | 52 ## @item |
10846
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
53 ## @var{m} = @var{m1} * @var{m2} is the (left) preconditioning matrix, so that |
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
54 ## the iteration is (theoretically) equivalent to solving by @code{pcg} |
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
55 ## @code{@var{P} * |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
56 ## @var{x} = @var{m} \ @var{b}}, with @code{@var{P} = @var{m} \ @var{A}}. |
5837 | 57 ## 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
|
58 ## improve the overall performance of the method. Instead of matrices |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
59 ## @var{m1} and @var{m2}, the user may pass two functions which return |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
60 ## the results of applying the inverse of @var{m1} and @var{m2} to |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
61 ## a vector (usually this is the preferred way of using the preconditioner). |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
62 ## If @code{[]} is supplied for @var{m1}, or @var{m1} is omitted, no |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
63 ## preconditioning is applied. If @var{m2} is omitted, @var{m} = @var{m1} |
6747 | 64 ## will be used as preconditioner. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
65 ## |
5837 | 66 ## @item |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
67 ## @var{x0} is the initial guess. If @var{x0} is empty or omitted, the |
5837 | 68 ## function sets @var{x0} to a zero vector by default. |
69 ## @end itemize | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
70 ## |
5837 | 71 ## 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:
10846
diff
changeset
|
72 ## 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
|
73 ## which are passed to @code{pcg}. See the examples below for further |
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
74 ## details. The output arguments are |
5837 | 75 ## |
76 ## @itemize | |
77 ## @item | |
78 ## @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:
10846
diff
changeset
|
79 ## @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
|
80 ## |
5837 | 81 ## @item |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
82 ## @var{flag} reports on the convergence. @code{@var{flag} = 0} means |
5837 | 83 ## 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
|
84 ## 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
|
85 ## for the iteration count was reached. @code{@var{flag} = 3} reports that |
5837 | 86 ## the (preconditioned) matrix was found not positive definite. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
87 ## |
5837 | 88 ## @item |
89 ## @var{relres} is the ratio of the final residual to its initial value, | |
90 ## 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
|
91 ## |
5837 | 92 ## @item |
93 ## @var{iter} is the actual number of iterations performed. | |
94 ## | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
95 ## @item |
5837 | 96 ## @var{resvec} describes the convergence history of the method. |
97 ## @code{@var{resvec} (i,1)} is the Euclidean norm of the residual, and | |
98 ## @code{@var{resvec} (i,2)} is the preconditioned residual norm, | |
99 ## after the (@var{i}-1)-th iteration, @code{@var{i} = | |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
100 ## 1, 2, @dots{}, @var{iter}+1}. The preconditioned residual norm |
6547 | 101 ## is defined as |
5838 | 102 ## @code{norm (@var{r}) ^ 2 = @var{r}' * (@var{m} \ @var{r})} where |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
103 ## @code{@var{r} = @var{b} - @var{A} * @var{x}}, see also the |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
104 ## description of @var{m}. If @var{eigest} is not required, only |
5837 | 105 ## @code{@var{resvec} (:,1)} is returned. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
106 ## |
5837 | 107 ## @item |
108 ## @var{eigest} returns the estimate for the smallest @code{@var{eigest} | |
109 ## (1)} and largest @code{@var{eigest} (2)} eigenvalues of the | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
110 ## preconditioned matrix @code{@var{P} = @var{m} \ @var{A}}. In |
7007 | 111 ## particular, if no preconditioning is used, the estimates for the |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
112 ## extreme eigenvalues of @var{A} are returned. @code{@var{eigest} (1)} |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
113 ## is an overestimate and @code{@var{eigest} (2)} is an underestimate, |
5837 | 114 ## so that @code{@var{eigest} (2) / @var{eigest} (1)} is a lower bound |
115 ## for @code{cond (@var{P}, 2)}, which nevertheless in the limit should | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
116 ## theoretically be equal to the actual value of the condition number. |
5837 | 117 ## The method which computes @var{eigest} works only for symmetric positive |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
118 ## definite @var{A} and @var{m}, and the user is responsible for |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
119 ## verifying this assumption. |
5837 | 120 ## @end itemize |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
121 ## |
5837 | 122 ## 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
|
123 ## sparsity of A) |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
124 ## |
5837 | 125 ## @example |
126 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
127 ## n = 10; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
128 ## A = diag (sparse (1:n)); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
129 ## b = rand (n, 1); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
130 ## [l, u, p, q] = luinc (A, 1.e-3); |
5837 | 131 ## @end group |
132 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
133 ## |
5837 | 134 ## @sc{Example 1:} Simplest use of @code{pcg} |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
135 ## |
5837 | 136 ## @example |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
137 ## x = pcg (A,b) |
5837 | 138 ## @end example |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
139 ## |
5837 | 140 ## @sc{Example 2:} @code{pcg} with a function which computes |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
141 ## @code{@var{A} * @var{x}} |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
142 ## |
5837 | 143 ## @example |
144 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
145 ## 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
|
146 ## y = [1:N]' .* x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
147 ## endfunction |
5837 | 148 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
149 ## x = pcg ("apply_a", b) |
5837 | 150 ## @end group |
151 ## @end example | |
6747 | 152 ## |
153 ## @sc{Example 3:} @code{pcg} with a preconditioner: @var{l} * @var{u} | |
154 ## | |
155 ## @example | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
156 ## x = pcg (A, b, 1.e-6, 500, l*u) |
6747 | 157 ## @end example |
158 ## | |
159 ## @sc{Example 4:} @code{pcg} with a preconditioner: @var{l} * @var{u}. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
160 ## Faster than @sc{Example 3} since lower and upper triangular matrices |
6747 | 161 ## are easier to invert |
162 ## | |
163 ## @example | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
164 ## x = pcg (A, b, 1.e-6, 500, l, u) |
6747 | 165 ## @end example |
166 ## | |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
167 ## @sc{Example 5:} Preconditioned iteration, with full diagnostics. The |
5837 | 168 ## 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:
10846
diff
changeset
|
169 ## @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
|
170 ## |
5837 | 171 ## @example |
172 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
173 ## 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
|
174 ## 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
|
175 ## y = x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
176 ## 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
|
177 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
178 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
179 ## [x, flag, relres, iter, resvec, eigest] = ... |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
180 ## pcg (A, b, [], [], "apply_m"); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
181 ## semilogy (1:iter+1, resvec); |
5837 | 182 ## @end group |
183 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
184 ## |
6747 | 185 ## @sc{Example 6:} Finally, a preconditioner which depends on a |
5837 | 186 ## parameter @var{k}. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
187 ## |
5837 | 188 ## @example |
189 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
190 ## function y = apply_M (x, varargin) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
191 ## K = varargin@{1@}; |
6547 | 192 ## y = x; |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
193 ## 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
|
194 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
195 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
196 ## [x, flag, relres, iter, resvec, eigest] = ... |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
197 ## pcg (A, b, [], [], "apply_m", [], [], 3) |
5837 | 198 ## @end group |
199 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
200 ## |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
201 ## References: |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
202 ## |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10791
diff
changeset
|
203 ## @enumerate |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10791
diff
changeset
|
204 ## @item |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10791
diff
changeset
|
205 ## C.T. Kelley, @cite{Iterative Methods for Linear and Nonlinear Equations}, |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10791
diff
changeset
|
206 ## SIAM, 1995. (the base PCG algorithm) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
207 ## |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10791
diff
changeset
|
208 ## @item |
16826
a4969508008e
doc: Periodic spellcheck of the documentation.
Rik <rik@octave.org>
parents:
14872
diff
changeset
|
209 ## Y. Saad, @cite{Iterative Methods for Sparse Linear Systems}, |
a4969508008e
doc: Periodic spellcheck of the documentation.
Rik <rik@octave.org>
parents:
14872
diff
changeset
|
210 ## @nospell{PWS} 1996. (condition number estimate from PCG) |
a4969508008e
doc: Periodic spellcheck of the documentation.
Rik <rik@octave.org>
parents:
14872
diff
changeset
|
211 ## Revised version of this book is available online at |
a4969508008e
doc: Periodic spellcheck of the documentation.
Rik <rik@octave.org>
parents:
14872
diff
changeset
|
212 ## @url{http://www-users.cs.umn.edu/~saad/books.html} |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10791
diff
changeset
|
213 ## @end enumerate |
5837 | 214 ## |
215 ## @seealso{sparse, pcr} | |
216 ## @end deftypefn | |
217 | |
5838 | 218 ## Author: Piotr Krzyzanowski <piotr.krzyzanowski@mimuw.edu.pl> |
6747 | 219 ## Modified by: Vittoria Rezzonico <vittoria.rezzonico@epfl.ch> |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
220 ## - Add the ability to provide the pre-conditioner as two separate matrices |
5837 | 221 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
222 function [x, flag, relres, iter, resvec, eigest] = pcg (A, b, tol, maxit, m1, m2, x0, varargin) |
5837 | 223 |
8507 | 224 ## M = M1*M2 |
6747 | 225 |
226 if (nargin < 7 || isempty (x0)) | |
5838 | 227 x = zeros (size (b)); |
5837 | 228 else |
229 x = x0; | |
230 endif | |
231 | |
8507 | 232 if (nargin < 5 || isempty (m1)) |
233 exist_m1 = 0; | |
234 else | |
235 exist_m1 = 1; | |
236 endif | |
6747 | 237 |
8507 | 238 if (nargin < 6 || isempty (m2)) |
239 exist_m2 = 0; | |
240 else | |
241 exist_m2 = 1; | |
242 endif | |
5837 | 243 |
5838 | 244 if (nargin < 4 || isempty (maxit)) |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14868
diff
changeset
|
245 maxit = min (rows (b), 20); |
5837 | 246 endif |
247 | |
5838 | 248 maxit += 2; |
5837 | 249 |
5838 | 250 if (nargin < 3 || isempty (tol)) |
5837 | 251 tol = 1e-6; |
252 endif | |
253 | |
254 preconditioned_residual_out = false; | |
255 if (nargout > 5) | |
5838 | 256 T = zeros (maxit, maxit); |
5837 | 257 preconditioned_residual_out = true; |
258 endif | |
259 | |
8506 | 260 ## Assume A is positive definite. |
261 matrix_positive_definite = true; | |
5837 | 262 |
5838 | 263 p = zeros (size (b)); |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
264 oldtau = 1; |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
265 if (isnumeric (A)) |
8506 | 266 ## A is a matrix. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
267 r = b - A*x; |
8506 | 268 else |
269 ## A should be a function. | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
270 r = b - feval (A, x, varargin{:}); |
5837 | 271 endif |
272 | |
5838 | 273 resvec(1,1) = norm (r); |
5837 | 274 alpha = 1; |
275 iter = 2; | |
276 | |
6747 | 277 while (resvec (iter-1,1) > tol * resvec (1,1) && iter < maxit) |
8507 | 278 if (exist_m1) |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
279 if (isnumeric (m1)) |
10549 | 280 y = m1 \ r; |
6747 | 281 else |
10549 | 282 y = feval (m1, r, varargin{:}); |
5837 | 283 endif |
6747 | 284 else |
285 y = r; | |
286 endif | |
8507 | 287 if (exist_m2) |
288 if (isnumeric (m2)) | |
10549 | 289 z = m2 \ y; |
6747 | 290 else |
10549 | 291 z = feval (m2, y, varargin{:}); |
6747 | 292 endif |
293 else | |
294 z = y; | |
5837 | 295 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
296 tau = z' * r; |
6747 | 297 resvec (iter-1,2) = sqrt (tau); |
5838 | 298 beta = tau / oldtau; |
5837 | 299 oldtau = tau; |
6747 | 300 p = z + beta * p; |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
301 if (isnumeric (A)) |
8506 | 302 ## A is a matrix. |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
303 w = A * p; |
8506 | 304 else |
305 ## A should be a function. | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10846
diff
changeset
|
306 w = feval (A, p, varargin{:}); |
5837 | 307 endif |
8506 | 308 ## Needed only for eigest. |
309 oldalpha = alpha; | |
5838 | 310 alpha = tau / (p'*w); |
8506 | 311 if (alpha <= 0.0) |
312 ## Negative matrix. | |
5837 | 313 matrix_positive_definite = false; |
314 endif | |
6747 | 315 x += alpha * p; |
316 r -= alpha * w; | |
5838 | 317 if (nargout > 5 && iter > 2) |
5837 | 318 T(iter-1:iter, iter-1:iter) = T(iter-1:iter, iter-1:iter) + ... |
10549 | 319 [1 sqrt(beta); sqrt(beta) beta]./oldalpha; |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
320 ## EVS = eig (T(2:iter-1,2:iter-1)); |
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
321 ## fprintf (stderr,"PCG condest: %g (iteration: %d)\n", max (EVS)/min (EVS),iter); |
5837 | 322 endif |
6747 | 323 resvec (iter,1) = norm (r); |
5838 | 324 iter++; |
5837 | 325 endwhile |
326 | |
327 if (nargout > 5) | |
5838 | 328 if (matrix_positive_definite) |
5837 | 329 if (iter > 3) |
10549 | 330 T = T(2:iter-2,2:iter-2); |
331 l = eig (T); | |
332 eigest = [min(l), max(l)]; | |
333 ## fprintf (stderr, "pcg condest: %g\n", eigest(2)/eigest(1)); | |
5837 | 334 else |
10549 | 335 eigest = [NaN, NaN]; |
11588
d5bd2766c640
style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
336 warning ("pcg: eigenvalue estimate failed: iteration converged too fast"); |
5837 | 337 endif |
338 else | |
5838 | 339 eigest = [NaN, NaN]; |
5837 | 340 endif |
341 | |
8506 | 342 ## Apply the preconditioner once more and finish with the precond |
343 ## residual. | |
8507 | 344 if (exist_m1) |
345 if (isnumeric (m1)) | |
10549 | 346 y = m1 \ r; |
6747 | 347 else |
10549 | 348 y = feval (m1, r, varargin{:}); |
5837 | 349 endif |
6747 | 350 else |
351 y = r; | |
5837 | 352 endif |
8507 | 353 if (exist_m2) |
354 if (isnumeric (m2)) | |
10549 | 355 z = m2 \ y; |
6747 | 356 else |
10549 | 357 z = feval (m2, y, varargin{:}); |
6747 | 358 endif |
359 else | |
360 z = y; | |
361 endif | |
362 | |
363 resvec (iter-1,2) = sqrt (r' * z); | |
5837 | 364 else |
5838 | 365 resvec = resvec(:,1); |
5837 | 366 endif |
367 | |
368 flag = 0; | |
6747 | 369 relres = resvec (iter-1,1) ./ resvec(1,1); |
5838 | 370 iter -= 2; |
6747 | 371 if (iter >= maxit - 2) |
5837 | 372 flag = 1; |
373 if (nargout < 2) | |
6498 | 374 warning ("pcg: maximum number of iterations (%d) reached\n", iter); |
6747 | 375 warning ("the initial residual norm was reduced %g times.\n", ... |
10549 | 376 1.0 / relres); |
5837 | 377 endif |
5838 | 378 elseif (nargout < 2) |
6498 | 379 fprintf (stderr, "pcg: converged in %d iterations. ", iter); |
380 fprintf (stderr, "the initial residual norm was reduced %g times.\n",... | |
10549 | 381 1.0/relres); |
5837 | 382 endif |
383 | |
5838 | 384 if (! matrix_positive_definite) |
5837 | 385 flag = 3; |
386 if (nargout < 2) | |
6498 | 387 warning ("pcg: matrix not positive definite?\n"); |
5837 | 388 endif |
389 endif | |
390 endfunction | |
391 | |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
392 |
5837 | 393 %!demo |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
394 %! ## Simplest usage of pcg (see also 'help pcg') |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
395 %! |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
396 %! N = 10; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
397 %! A = diag ([1:N]); b = rand (N, 1); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
398 %! y = A \ b; # y is the true solution |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
399 %! x = pcg (A, b); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
400 %! printf ("The solution relative error is %g\n", norm (x - y) / norm (y)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
401 %! |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
402 %! ## You shouldn't be afraid if pcg issues some warning messages in this |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
403 %! ## example: watch out in the second example, why it takes N iterations |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
404 %! ## of pcg to converge to (a very accurate, by the way) solution |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
405 |
5837 | 406 %!demo |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
407 %! ## Full output from pcg, except for the eigenvalue estimates |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
408 %! ## We use this output to plot the convergence history |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
409 %! |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
410 %! N = 10; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
411 %! A = diag ([1:N]); b = rand (N, 1); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
412 %! X = A \ b; # X is the true solution |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
413 %! [x, flag, relres, iter, resvec] = pcg (A, b); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
414 %! printf ("The solution relative error is %g\n", norm (x - X) / norm (X)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
415 %! title ("Convergence history"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
416 %! semilogy ([0:iter], resvec / resvec(1), "o-g"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
417 %! xlabel ("Iteration"); ylabel ("log(||b-Ax||/||b||)"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
418 %! legend ("relative residual"); |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
419 |
5837 | 420 %!demo |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
421 %! ## Full output from pcg, including the eigenvalue estimates |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
422 %! ## Hilbert matrix is extremely ill-conditioned, so pcg WILL have problems |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
423 %! |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
424 %! N = 10; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
425 %! A = hilb (N); b = rand (N, 1); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
426 %! X = A \ b; # X is the true solution |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
427 %! [x, flag, relres, iter, resvec, eigest] = pcg (A, b, [], 200); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
428 %! printf ("The solution relative error is %g\n", norm (x - X) / norm (X)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
429 %! printf ("Condition number estimate is %g\n", eigest(2) / eigest(1)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
430 %! printf ("Actual condition number is %g\n", cond (A)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
431 %! title ("Convergence history"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
432 %! semilogy ([0:iter], resvec, ["o-g";"+-r"]); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
433 %! xlabel ("Iteration"); ylabel ("log(||b-Ax||)"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
434 %! legend ("absolute residual", "absolute preconditioned residual"); |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
435 |
5837 | 436 %!demo |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
437 %! ## Full output from pcg, including the eigenvalue estimates |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
438 %! ## We use the 1-D Laplacian matrix for A, and cond(A) = O(N^2) |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
439 %! ## and that's the reason we need some preconditioner; here we take |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
440 %! ## a very simple and not powerful Jacobi preconditioner, |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
441 %! ## which is the diagonal of A. |
5837 | 442 %! |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
443 %! N = 100; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
444 %! A = zeros (N, N); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
445 %! for i = 1 : N - 1 # form 1-D Laplacian matrix |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
446 %! A(i:i+1, i:i+1) = [2 -1; -1 2]; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
447 %! endfor |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
448 %! b = rand (N, 1); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
449 %! X = A \ b; # X is the true solution |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
450 %! maxit = 80; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
451 %! printf ("System condition number is %g\n", cond (A)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
452 %! ## No preconditioner: the convergence is very slow! |
5837 | 453 %! |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
454 %! [x, flag, relres, iter, resvec, eigest] = pcg (A, b, [], maxit); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
455 %! printf ("System condition number estimate is %g\n", eigest(2) / eigest(1)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
456 %! title ("Convergence history"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
457 %! semilogy ([0:iter], resvec(:,1), "o-g"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
458 %! xlabel ("Iteration"); ylabel ("log(||b-Ax||)"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
459 %! legend ("NO preconditioning: absolute residual"); |
5837 | 460 %! |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
461 %! pause (1); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
462 %! ## Test Jacobi preconditioner: it will not help much!!! |
5837 | 463 %! |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
464 %! M = diag (diag (A)); # Jacobi preconditioner |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
465 %! [x, flag, relres, iter, resvec, eigest] = pcg (A, b, [], maxit, M); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
466 %! printf ("JACOBI preconditioned system condition number estimate is %g\n", eigest(2) / eigest(1)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
467 %! hold on; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
468 %! semilogy ([0:iter], resvec(:,1), "o-r"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
469 %! legend ("NO preconditioning: absolute residual", ... |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
470 %! "JACOBI preconditioner: absolute residual"); |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
471 %! |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
472 %! pause (1); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
473 %! ## Test nonoverlapping block Jacobi preconditioner: it will help much! |
5837 | 474 %! |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
475 %! M = zeros (N, N); k = 4; |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
476 %! for i = 1 : k : N # form 1-D Laplacian matrix |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
477 %! M(i:i+k-1, i:i+k-1) = A(i:i+k-1, i:i+k-1); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
478 %! endfor |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
479 %! [x, flag, relres, iter, resvec, eigest] = pcg (A, b, [], maxit, M); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
480 %! printf ("BLOCK JACOBI preconditioned system condition number estimate is %g\n", eigest(2) / eigest(1)); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
481 %! semilogy ([0:iter], resvec(:,1), "o-b"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
482 %! legend ("NO preconditioning: absolute residual", ... |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
483 %! "JACOBI preconditioner: absolute residual", ... |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
484 %! "BLOCK JACOBI preconditioner: absolute residual"); |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
485 %! hold off; |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
486 |
5837 | 487 %!test |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
488 %! ## solve small diagonal system |
5837 | 489 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
490 %! N = 10; |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
491 %! A = diag ([1:N]); b = rand (N, 1); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
492 %! X = A \ b; # X is the true solution |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
493 %! [x, flag] = pcg (A, b, [], N+1); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
494 %! assert (norm (x - X) / norm (X), 0, 1e-10); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
495 %! assert (flag, 0); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
496 |
5837 | 497 %!test |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
498 %! ## solve small indefinite diagonal system |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
499 %! ## despite A is indefinite, the iteration continues and converges |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
500 %! ## indefiniteness of A is detected |
5837 | 501 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
502 %! N = 10; |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
503 %! A = diag([1:N] .* (-ones(1, N) .^ 2)); b = rand (N, 1); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
504 %! X = A \ b; # X is the true solution |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
505 %! [x, flag] = pcg (A, b, [], N+1); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
506 %! assert (norm (x - X) / norm (X), 0, 1e-10); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
507 %! assert (flag, 3); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
508 |
5837 | 509 %!test |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
510 %! ## solve tridiagonal system, do not converge in default 20 iterations |
5837 | 511 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
512 %! N = 100; |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
513 %! A = zeros (N, N); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
514 %! for i = 1 : N - 1 # form 1-D Laplacian matrix |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
515 %! A(i:i+1, i:i+1) = [2 -1; -1 2]; |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
516 %! endfor |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
517 %! b = ones (N, 1); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
518 %! X = A \ b; # X is the true solution |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
519 %! [x, flag, relres, iter, resvec, eigest] = pcg (A, b, 1e-12); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
520 %! assert (flag); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
521 %! assert (relres > 1.0); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
522 %! assert (iter, 20); # should perform max allowable default number of iterations |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
523 |
5837 | 524 %!test |
17344
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
525 %! ## solve tridiagonal system with 'perfect' preconditioner |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
526 %! ## which converges in one iteration, so the eigest does not |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
527 %! ## work and issues a warning |
5837 | 528 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
529 %! N = 100; |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
530 %! A = zeros (N, N); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
531 %! for i = 1 : N - 1 # form 1-D Laplacian matrix |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
532 %! A (i:i+1, i:i+1) = [2 -1; -1 2]; |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
533 %! endfor |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
534 %! b = ones (N, 1); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
535 %! X = A \ b; # X is the true solution |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
536 %! [x, flag, relres, iter, resvec, eigest] = pcg (A, b, [], [], A, [], b); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
537 %! assert (norm (x - X) / norm (X), 0, 1e-6); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
538 %! assert (flag, 0); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
539 %! assert (iter, 1); # should converge in one iteration |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
540 %! assert (isnan (eigest), isnan ([NaN, NaN])); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
541 |