Mercurial > octave
annotate scripts/sparse/pcr.m @ 27919:1891570abac8
update Octave Project Developers copyright for the new year
In files that have the "Octave Project Developers" copyright notice,
update for 2020.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Mon, 06 Jan 2020 22:29:51 -0500 |
parents | b442ec6dda5c |
children | bd51beb6205e |
rev | line source |
---|---|
27919
1891570abac8
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
27918
diff
changeset
|
1 ## Copyright (C) 2004-2020 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
2 ## |
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
3 ## See the file COPYRIGHT.md in the top-level directory of this distribution |
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
4 ## or <https://octave.org/COPYRIGHT.html/>. |
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
5 ## |
5837 | 6 ## |
7 ## This file is part of Octave. | |
8 ## | |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
9 ## Octave is free software: you can redistribute it and/or modify it |
5837 | 10 ## under the terms of the GNU General Public License as published by |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
11 ## the Free Software Foundation, either version 3 of the License, or |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
12 ## (at your option) any later version. |
5837 | 13 ## |
14 ## Octave is distributed in the hope that it will be useful, but | |
15 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
16 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
17 ## GNU General Public License for more details. |
5837 | 18 ## |
19 ## You should have received a copy of the GNU General Public License | |
7016 | 20 ## along with Octave; see the file COPYING. If not, see |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
21 ## <https://www.gnu.org/licenses/>. |
5837 | 22 |
23 ## -*- texinfo -*- | |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
24 ## @deftypefn {} {@var{x} =} pcr (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{m}, @var{x0}, @dots{}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
25 ## @deftypefnx {} {[@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
|
26 ## |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
27 ## Solve the linear system of equations @code{@var{A} * @var{x} = @var{b}} by |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
28 ## means of the Preconditioned Conjugate Residuals iterative method. |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
29 ## |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
30 ## The input arguments are |
5837 | 31 ## |
32 ## @itemize | |
33 ## @item | |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
34 ## @var{A} can be either a square (preferably sparse) matrix or a function |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
35 ## handle, inline function or string containing the name of a function which |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
36 ## computes @code{@var{A} * @var{x}}. In principle @var{A} should be |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
37 ## symmetric and non-singular; if @code{pcr} finds @var{A} to be numerically |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
38 ## singular, you will get a warning message and the @var{flag} output |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
39 ## parameter will be set. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
40 ## |
5837 | 41 ## @item |
42 ## @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
|
43 ## |
5837 | 44 ## @item |
45 ## @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
|
46 ## @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
|
47 ## @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
|
48 ## @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
|
49 ## If @var{tol} is empty or is omitted, the function sets |
5837 | 50 ## @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
|
51 ## |
5837 | 52 ## @item |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
53 ## @var{maxit} is the maximum allowable number of iterations; if @code{[]} is |
25472
9771111f04f4
doc: Use @var rather than @code to mark inputs to functions in docstrings.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
54 ## supplied for @var{maxit}, or @code{pcr} has less arguments, a default |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
55 ## value equal to 20 is used. |
5837 | 56 ## |
57 ## @item | |
5838 | 58 ## @var{m} is the (left) preconditioning matrix, so that the iteration is |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
59 ## (theoretically) equivalent to solving by |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
60 ## @code{pcr} @code{@var{P} * @var{x} = @var{m} \ @var{b}}, with |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
61 ## @code{@var{P} = @var{m} \ @var{A}}. Note that a proper choice of the |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
62 ## preconditioner may dramatically improve the overall performance of the |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
63 ## method. Instead of matrix @var{m}, the user may pass a function which |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
64 ## returns the results of applying the inverse of @var{m} to a vector |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
65 ## (usually this is the preferred way of using the preconditioner). If |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
66 ## @code{[]} is supplied for @var{m}, or @var{m} is omitted, no |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
67 ## preconditioning is applied. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
68 ## |
5837 | 69 ## @item |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
70 ## @var{x0} is the initial guess. If @var{x0} is empty or omitted, the |
5837 | 71 ## function sets @var{x0} to a zero vector by default. |
72 ## @end itemize | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
73 ## |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
74 ## The arguments which follow @var{x0} are treated as parameters, and passed |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
75 ## in a proper way to any of the functions (@var{A} or @var{m}) which are |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
76 ## passed to @code{pcr}. See the examples below for further details. |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
77 ## |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
78 ## The output arguments are |
5837 | 79 ## |
80 ## @itemize | |
81 ## @item | |
82 ## @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
|
83 ## @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
|
84 ## |
5837 | 85 ## @item |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
86 ## @var{flag} reports on the convergence. @code{@var{flag} = 0} means the |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
87 ## solution converged and the tolerance criterion given by @var{tol} is |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
88 ## satisfied. @code{@var{flag} = 1} means that the @var{maxit} limit for the |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
89 ## iteration count was reached. @code{@var{flag} = 3} reports a @code{pcr} |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
90 ## 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
|
91 ## |
5837 | 92 ## @item |
93 ## @var{relres} is the ratio of the final residual to its initial value, | |
94 ## 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
|
95 ## |
5837 | 96 ## @item |
97 ## @var{iter} is the actual number of iterations performed. | |
98 ## | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
99 ## @item |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
100 ## @var{resvec} describes the convergence history of the method, so that |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
101 ## @code{@var{resvec} (i)} contains the Euclidean norms of the residual after |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
102 ## the (@var{i}-1)-th iteration, @code{@var{i} = 1,2, @dots{}, @var{iter}+1}. |
5837 | 103 ## @end itemize |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
104 ## |
5837 | 105 ## 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
|
106 ## sparsity of A) |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
107 ## |
5837 | 108 ## @example |
109 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
110 ## n = 10; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
111 ## 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
|
112 ## b = rand (N, 1); |
5837 | 113 ## @end group |
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 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
|
117 ## |
5837 | 118 ## @example |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
119 ## x = pcr (A, b) |
5837 | 120 ## @end example |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
121 ## |
5837 | 122 ## @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
|
123 ## @code{@var{A} * @var{x}}. |
5837 | 124 ## |
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 ## 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
|
128 ## y = [1:10]' .* x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
129 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
130 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
131 ## x = pcr ("apply_a", b) |
5837 | 132 ## @end group |
133 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
134 ## |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
135 ## @sc{Example 3:} Preconditioned iteration, with full diagnostics. The |
5837 | 136 ## 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
|
137 ## @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
|
138 ## |
5837 | 139 ## @example |
140 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
141 ## 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
|
142 ## 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
|
143 ## y = x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
144 ## 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
|
145 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
146 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
147 ## [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
|
148 ## 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
|
149 ## semilogy ([1:iter+1], resvec); |
5837 | 150 ## @end group |
151 ## @end example | |
152 ## | |
153 ## @sc{Example 4:} Finally, a preconditioner which depends on a | |
154 ## parameter @var{k}. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
155 ## |
5837 | 156 ## @example |
157 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
158 ## 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
|
159 ## k = varargin@{1@}; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
160 ## y = x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
161 ## 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
|
162 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
163 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
164 ## [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
|
165 ## pcr (A, b, [], [], "apply_m"', [], 3) |
5837 | 166 ## @end group |
167 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
168 ## |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
169 ## References: |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
170 ## |
19040
0850b5212619
doc: Add @nospell macro around proper names in documentation.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
171 ## [1] @nospell{W. Hackbusch}, @cite{Iterative Solution of Large Sparse |
21315
ea2c2e08ceda
doc: Periodic spellcheck of documentation.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
172 ## Systems of Equations}, section 9.5.4; @nospell{Springer}, 1994 |
5837 | 173 ## |
174 ## @seealso{sparse, pcg} | |
175 ## @end deftypefn | |
176 | |
5838 | 177 ## Author: Piotr Krzyzanowski <piotr.krzyzanowski@mimuw.edu.pl> |
5837 | 178 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
179 function [x, flag, relres, iter, resvec] = pcr (A, b, tol, maxit, m, x0, varargin) |
5837 | 180 |
181 breakdown = false; | |
182 | |
5838 | 183 if (nargin < 6 || isempty (x0)) |
184 x = zeros (size (b)); | |
5837 | 185 else |
186 x = x0; | |
187 endif | |
188 | |
189 if (nargin < 5) | |
8507 | 190 m = []; |
5837 | 191 endif |
192 | |
5838 | 193 if (nargin < 4 || isempty (maxit)) |
5837 | 194 maxit = 20; |
195 endif | |
196 | |
5838 | 197 maxit += 2; |
5837 | 198 |
5838 | 199 if (nargin < 3 || isempty (tol)) |
5837 | 200 tol = 1e-6; |
201 endif | |
202 | |
203 if (nargin < 2) | |
5838 | 204 print_usage (); |
5837 | 205 endif |
206 | |
207 ## init | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
208 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
|
209 r = b - A*x; |
10549 | 210 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
|
211 r = b - feval (A, x, varargin{:}); |
5837 | 212 endif |
213 | |
10549 | 214 if (isnumeric (m)) # is M a matrix? |
215 if (isempty (m)) # if M is empty, use no precond | |
5837 | 216 p = r; |
10549 | 217 else # otherwise, apply the precond |
8507 | 218 p = m \ r; |
5837 | 219 endif |
10549 | 220 else # then M should be a function! |
8507 | 221 p = feval (m, r, varargin{:}); |
5837 | 222 endif |
223 | |
224 iter = 2; | |
225 | |
226 b_bot_old = 1; | |
5838 | 227 q_old = p_old = s_old = zeros (size (x)); |
5837 | 228 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
229 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
|
230 q = A * p; |
10549 | 231 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
|
232 q = feval (A, p, varargin{:}); |
5837 | 233 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
234 |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
235 resvec(1) = abs (norm (r)); |
5837 | 236 |
237 ## iteration | |
5838 | 238 while (resvec(iter-1) > tol*resvec(1) && iter < maxit) |
239 | |
10549 | 240 if (isnumeric (m)) # is M a matrix? |
241 if (isempty (m)) # if M is empty, use no precond | |
242 s = q; | |
243 else # otherwise, apply the precond | |
244 s = m \ q; | |
5837 | 245 endif |
10549 | 246 else # then M should be a function! |
8507 | 247 s = feval (m, q, varargin{:}); |
5837 | 248 endif |
5838 | 249 b_top = r' * s; |
250 b_bot = q' * s; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
251 |
5837 | 252 if (b_bot == 0.0) |
253 breakdown = true; | |
254 break; | |
255 endif | |
5838 | 256 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
|
257 |
5838 | 258 x += lambda*p; |
259 r -= lambda*q; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
260 |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
261 if (isnumeric (A)) # is A a matrix? |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
262 t = A*s; |
10549 | 263 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
|
264 t = feval (A, s, varargin{:}); |
5837 | 265 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
266 |
5838 | 267 alpha0 = (t'*s) / b_bot; |
268 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
|
269 |
5838 | 270 p_temp = p; |
271 q_temp = q; | |
272 | |
5837 | 273 p = s - alpha0*p - alpha1*p_old; |
274 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
|
275 |
5837 | 276 s_old = s; |
277 p_old = p_temp; | |
278 q_old = q_temp; | |
279 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
|
280 |
5838 | 281 resvec(iter) = abs (norm (r)); |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20727
diff
changeset
|
282 iter += 1; |
5837 | 283 endwhile |
284 | |
285 flag = 0; | |
5838 | 286 relres = resvec(iter-1) ./ resvec(1); |
287 iter -= 2; | |
288 if (iter >= maxit-2) | |
5837 | 289 flag = 1; |
290 if (nargout < 2) | |
6498 | 291 warning ("pcr: maximum number of iterations (%d) reached\n", iter); |
20727
a5949b3d2332
Preface warning() messages with name of function when possible.
Rik <rik@octave.org>
parents:
20164
diff
changeset
|
292 warning ("pcr: the initial residual norm was reduced %g times\n", |
a5949b3d2332
Preface warning() messages with name of function when possible.
Rik <rik@octave.org>
parents:
20164
diff
changeset
|
293 1.0/relres); |
5837 | 294 endif |
5838 | 295 elseif (nargout < 2 && ! breakdown) |
6498 | 296 fprintf (stderr, "pcr: converged in %d iterations. \n", iter); |
20727
a5949b3d2332
Preface warning() messages with name of function when possible.
Rik <rik@octave.org>
parents:
20164
diff
changeset
|
297 fprintf (stderr, "pcr: the initial residual norm was reduced %g times\n", |
10549 | 298 1.0 / relres); |
5837 | 299 endif |
5838 | 300 |
5837 | 301 if (breakdown) |
302 flag = 3; | |
303 if (nargout < 2) | |
7001 | 304 warning ("pcr: breakdown occurred:\n"); |
6498 | 305 warning ("system matrix singular or preconditioner indefinite?\n"); |
5837 | 306 endif |
307 endif | |
308 | |
309 endfunction | |
310 | |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
311 |
5837 | 312 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
313 %! ## Simplest usage of PCR (see also 'help pcr') |
5837 | 314 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
315 %! N = 20; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
316 %! A = diag (linspace (-3.1,3,N)); b = rand (N,1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
317 %! y = A \ b; # y is the true solution |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
318 %! x = pcr (A,b); |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
319 %! printf ("The solution relative error is %g\n", norm (x-y) / norm (y)); |
5837 | 320 %! |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
321 %! ## You shouldn't be afraid if PCR issues some warning messages in this |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
322 %! ## 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:
14868
diff
changeset
|
323 %! ## of PCR 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
|
324 |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
325 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
326 %! ## Full output from PCR |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
327 %! ## We use this output to plot the convergence history |
5837 | 328 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
329 %! N = 20; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
330 %! A = diag (linspace (-3.1,30,N)); b = rand (N,1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
331 %! X = A \ b; # X is the true solution |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
332 %! [x, flag, relres, iter, resvec] = pcr (A,b); |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
333 %! printf ("The solution relative error is %g\n", norm (x-X) / norm (X)); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
334 %! clf; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
335 %! title ("Convergence history"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
336 %! xlabel ("Iteration"); ylabel ("log(||b-Ax||/||b||)"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
337 %! semilogy ([0:iter], resvec/resvec(1), "o-g;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
|
338 |
5837 | 339 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
340 %! ## Full output from PCR |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
341 %! ## We use indefinite matrix based on the Hilbert matrix, with one |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
342 %! ## strongly negative eigenvalue |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
343 %! ## Hilbert matrix is extremely ill conditioned, so is ours, |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
344 %! ## and that's why PCR WILL have problems |
5837 | 345 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
346 %! N = 10; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
347 %! A = hilb (N); A(1,1) = -A(1,1); b = rand (N,1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
348 %! X = A \ b; # X is the true solution |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
349 %! printf ("Condition number of A is %g\n", cond (A)); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
350 %! [x, flag, relres, iter, resvec] = pcr (A,b,[],200); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
351 %! if (flag == 3) |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
352 %! printf ("PCR breakdown. System matrix is [close to] singular\n"); |
22302
1c4cd12987f5
Use Octave syntax in graphics demos.
Rik <rik@octave.org>
parents:
21634
diff
changeset
|
353 %! endif |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
354 %! clf; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
355 %! title ("Convergence history"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
356 %! xlabel ("Iteration"); ylabel ("log(||b-Ax||)"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
357 %! semilogy ([0:iter], resvec, "o-g;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
|
358 |
5837 | 359 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
360 %! ## Full output from PCR |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
361 %! ## We use an indefinite matrix based on the 1-D Laplacian matrix for A, |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
362 %! ## and here we have cond(A) = O(N^2) |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
363 %! ## 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:
14868
diff
changeset
|
364 %! ## a very simple and not powerful Jacobi preconditioner, |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
365 %! ## which is the diagonal of A. |
5837 | 366 %! |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
367 %! ## Note that we use here indefinite preconditioners! |
5837 | 368 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
369 %! N = 100; |
25755
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
370 %! ## Form 1-D Laplacian matrix |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
371 %! A = 2 * eye (N,N); |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
372 %! A(2:(N+1):end) = -1; |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
373 %! A((N+1):(N+1):end) = -1; |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
374 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
375 %! A = [A, zeros(size(A)); zeros(size(A)), -A]; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
376 %! b = rand (2*N,1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
377 %! X = A \ b; # X is the true solution |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
378 %! maxit = 80; |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
379 %! printf ("System condition number is %g\n", cond (A)); |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
380 %! ## No preconditioner: the convergence is very slow! |
5837 | 381 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
382 %! [x, flag, relres, iter, resvec] = pcr (A,b,[],maxit); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
383 %! clf; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
384 %! title ("Convergence history"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
385 %! xlabel ("Iteration"); ylabel ("log(||b-Ax||)"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
386 %! semilogy ([0:iter], resvec, "o-g;NO preconditioning: absolute residual;"); |
5837 | 387 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
388 %! pause (1); |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
389 %! ## Test Jacobi preconditioner: it will not help much!!! |
5837 | 390 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
391 %! M = diag (diag (A)); # Jacobi preconditioner |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
392 %! [x, flag, relres, iter, resvec] = pcr (A,b,[],maxit,M); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
393 %! hold on; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
394 %! semilogy ([0:iter],resvec,"o-r;JACOBI preconditioner: absolute residual;"); |
5837 | 395 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
396 %! pause (1); |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
397 %! ## Test nonoverlapping block Jacobi preconditioner: this one should give |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
398 %! ## some convergence speedup! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
399 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
400 %! M = zeros (N,N); k = 4; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
401 %! for i=1:k:N # get k x k diagonal blocks of A |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
402 %! M(i:i+k-1,i:i+k-1) = A(i:i+k-1,i:i+k-1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
403 %! endfor |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
404 %! M = [M, zeros(size (M)); zeros(size(M)), -M]; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
405 %! [x, flag, relres, iter, resvec] = pcr (A,b,[],maxit,M); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
406 %! semilogy ([0:iter], resvec, "o-b;BLOCK JACOBI preconditioner: absolute residual;"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
407 %! 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
|
408 |
5837 | 409 %!test |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
410 %! ## solve small indefinite diagonal system |
5837 | 411 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
412 %! 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
|
413 %! A = diag (linspace (-10.1,10,N)); 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
|
414 %! 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
|
415 %! [x, flag] = pcr (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
|
416 %! assert (norm (x-X) / norm (X) < 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
|
417 %! 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
|
418 |
5837 | 419 %!test |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
420 %! ## solve tridiagonal system, do not converge in default 20 iterations |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
421 %! ## should perform max allowable default number of iterations |
5837 | 422 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
423 %! N = 100; |
25755
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
424 %! ## Form 1-D Laplacian matrix |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
425 %! A = 2 * eye (N,N); |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
426 %! A(2:(N+1):end) = -1; |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
427 %! A((N+1):(N+1):end) = -1; |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
428 %! 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
|
429 %! 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
|
430 %! [x, flag, relres, iter, resvec] = pcr (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
|
431 %! assert (flag, 1); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
432 %! assert (relres > 0.6); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
433 %! assert (iter, 20); |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
434 |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
435 %!test |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
436 %! ## solve tridiagonal system with "perfect" preconditioner |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
437 %! ## converges in one iteration |
5837 | 438 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
439 %! N = 100; |
25755
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
440 %! ## Form 1-D Laplacian matrix |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
441 %! A = 2 * eye (N,N); |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
442 %! A(2:(N+1):end) = -1; |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
443 %! A((N+1):(N+1):end) = -1; |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
444 %! 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
|
445 %! 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
|
446 %! [x, flag, relres, iter] = pcr (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
|
447 %! assert (norm (x-X) / norm(X) < 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
|
448 %! assert (relres < 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
|
449 %! 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
|
450 %! assert (iter, 1); # should converge in one iteration |