Mercurial > octave
annotate scripts/sparse/pcr.m @ 27984:b09432b20a84
maint: Remove special cases of old version control keywords in code base.
* common_size.m, cumtrapz.m, interp1.m, polyarea.m, trapz.m, inpolygon.m,
summer.m, javamem.m, duplication_matrix.m, subspace.m, starting_stepsize.m,
text.m, errorbar.m, loglogerr.m, peaks.m, __errplot__.m, semilogxerr.m,
semilogyerr.m, stemleaf.m, streamtube.m, __pltopt__.m, printd.m, pchip.m,
arch_fit.m, arch_rnd.m, arch_test.m, arma_rnd.m, autoreg_matrix.m, bartlett.m,
blackman.m, diffpara.m, durbinlevinson.m, fractdiff.m, hamming.m, hanning.m,
hurst.m, ifftshift.m, periodogram.m, rectangle_lw.m, rectangle_sw.m,
triangle_lw.m, triangle_sw.m, sinetone.m, sinewave.m, spectral_adf.m,
spectral_xdf.m, spencer.m, stft.m, synthesis.m, yulewalker.m, bicg.m,
bicgstab.m, cgs.m, gmres.m, pcr.m, sprand.m, tfqmr.m, betaincinv.m, betaln.m,
gammaincinv.m, hadamard.m, center.m, cov.m, discrete_inv.m, discrete_pdf.m,
discrete_rnd.m, empirical_cdf.m, empirical_inv.m, empirical_pdf.m,
empirical_rnd.m, iqr.m, kendall.m, mean.m, meansq.m, moment.m, prctile.m,
quantile.m, range.m, run_count.m, spearman.m, statistics.m, var.m, zscore.m,
datenum.m, datevec.m:
Remove special cases of old version control keywords in code base.
author | Rik <rik@octave.org> |
---|---|
date | Tue, 21 Jan 2020 14:04:31 -0800 |
parents | a4268efb7334 |
children | 0de38a6ef693 0a5b15007766 |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 ######################################################################## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 ## |
27919
1891570abac8
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
27918
diff
changeset
|
3 ## 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
|
4 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 ## See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 ## distribution or <https://octave.org/copyright/>. |
5837 | 7 ## |
8 ## This file is part of Octave. | |
9 ## | |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
10 ## Octave is free software: you can redistribute it and/or modify it |
5837 | 11 ## 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
|
12 ## 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
|
13 ## (at your option) any later version. |
5837 | 14 ## |
15 ## Octave is distributed in the hope that it will be useful, but | |
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
18 ## GNU General Public License for more details. |
5837 | 19 ## |
20 ## You should have received a copy of the GNU General Public License | |
7016 | 21 ## 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
|
22 ## <https://www.gnu.org/licenses/>. |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
23 ## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 ######################################################################## |
5837 | 25 |
26 ## -*- texinfo -*- | |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
27 ## @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
|
28 ## @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
|
29 ## |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
30 ## 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
|
31 ## 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
|
32 ## |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
33 ## The input arguments are |
5837 | 34 ## |
35 ## @itemize | |
36 ## @item | |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
37 ## @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
|
38 ## 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
|
39 ## 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
|
40 ## 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
|
41 ## 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
|
42 ## parameter will be set. |
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{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
|
46 ## |
5837 | 47 ## @item |
48 ## @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
|
49 ## @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
|
50 ## @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
|
51 ## @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
|
52 ## If @var{tol} is empty or is omitted, the function sets |
5837 | 53 ## @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
|
54 ## |
5837 | 55 ## @item |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
56 ## @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
|
57 ## 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
|
58 ## value equal to 20 is used. |
5837 | 59 ## |
60 ## @item | |
5838 | 61 ## @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
|
62 ## (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
|
63 ## @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
|
64 ## @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
|
65 ## 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
|
66 ## 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
|
67 ## 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
|
68 ## (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
|
69 ## @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
|
70 ## preconditioning is applied. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
71 ## |
5837 | 72 ## @item |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
73 ## @var{x0} is the initial guess. If @var{x0} is empty or omitted, the |
5837 | 74 ## function sets @var{x0} to a zero vector by default. |
75 ## @end itemize | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
76 ## |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
77 ## 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
|
78 ## 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
|
79 ## 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
|
80 ## |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
81 ## The output arguments are |
5837 | 82 ## |
83 ## @itemize | |
84 ## @item | |
85 ## @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
|
86 ## @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
|
87 ## |
5837 | 88 ## @item |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
89 ## @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
|
90 ## 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
|
91 ## 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
|
92 ## 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
|
93 ## 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
|
94 ## |
5837 | 95 ## @item |
96 ## @var{relres} is the ratio of the final residual to its initial value, | |
97 ## 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
|
98 ## |
5837 | 99 ## @item |
100 ## @var{iter} is the actual number of iterations performed. | |
101 ## | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
102 ## @item |
20164
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
103 ## @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
|
104 ## @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
|
105 ## the (@var{i}-1)-th iteration, @code{@var{i} = 1,2, @dots{}, @var{iter}+1}. |
5837 | 106 ## @end itemize |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
107 ## |
5837 | 108 ## 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
|
109 ## sparsity of A) |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
110 ## |
5837 | 111 ## @example |
112 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
113 ## n = 10; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
114 ## 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
|
115 ## b = rand (N, 1); |
5837 | 116 ## @end group |
117 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
118 ## |
5837 | 119 ## @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
|
120 ## |
5837 | 121 ## @example |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
122 ## x = pcr (A, b) |
5837 | 123 ## @end example |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
124 ## |
5837 | 125 ## @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
|
126 ## @code{@var{A} * @var{x}}. |
5837 | 127 ## |
128 ## @example | |
129 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
130 ## 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
|
131 ## y = [1:10]' .* x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
132 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
133 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
134 ## x = pcr ("apply_a", b) |
5837 | 135 ## @end group |
136 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
137 ## |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
138 ## @sc{Example 3:} Preconditioned iteration, with full diagnostics. The |
5837 | 139 ## 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
|
140 ## @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
|
141 ## |
5837 | 142 ## @example |
143 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
144 ## 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
|
145 ## 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
|
146 ## y = x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
147 ## 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
|
148 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
149 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
150 ## [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
|
151 ## 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
|
152 ## semilogy ([1:iter+1], resvec); |
5837 | 153 ## @end group |
154 ## @end example | |
155 ## | |
156 ## @sc{Example 4:} Finally, a preconditioner which depends on a | |
157 ## parameter @var{k}. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
158 ## |
5837 | 159 ## @example |
160 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
161 ## 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
|
162 ## k = varargin@{1@}; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
163 ## y = x; |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
164 ## 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
|
165 ## endfunction |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
166 ## |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
167 ## [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
|
168 ## pcr (A, b, [], [], "apply_m"', [], 3) |
5837 | 169 ## @end group |
170 ## @end example | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
171 ## |
27984
b09432b20a84
maint: Remove special cases of old version control keywords in code base.
Rik <rik@octave.org>
parents:
27978
diff
changeset
|
172 ## Reference: |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
173 ## |
27984
b09432b20a84
maint: Remove special cases of old version control keywords in code base.
Rik <rik@octave.org>
parents:
27978
diff
changeset
|
174 ## @nospell{W. Hackbusch}, @cite{Iterative Solution of Large Sparse |
21315
ea2c2e08ceda
doc: Periodic spellcheck of documentation.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
175 ## Systems of Equations}, section 9.5.4; @nospell{Springer}, 1994 |
5837 | 176 ## |
177 ## @seealso{sparse, pcg} | |
178 ## @end deftypefn | |
179 | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
180 function [x, flag, relres, iter, resvec] = pcr (A, b, tol, maxit, m, x0, varargin) |
5837 | 181 |
182 breakdown = false; | |
183 | |
5838 | 184 if (nargin < 6 || isempty (x0)) |
185 x = zeros (size (b)); | |
5837 | 186 else |
187 x = x0; | |
188 endif | |
189 | |
190 if (nargin < 5) | |
8507 | 191 m = []; |
5837 | 192 endif |
193 | |
5838 | 194 if (nargin < 4 || isempty (maxit)) |
5837 | 195 maxit = 20; |
196 endif | |
197 | |
5838 | 198 maxit += 2; |
5837 | 199 |
5838 | 200 if (nargin < 3 || isempty (tol)) |
5837 | 201 tol = 1e-6; |
202 endif | |
203 | |
204 if (nargin < 2) | |
5838 | 205 print_usage (); |
5837 | 206 endif |
207 | |
208 ## init | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
209 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
|
210 r = b - A*x; |
10549 | 211 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
|
212 r = b - feval (A, x, varargin{:}); |
5837 | 213 endif |
214 | |
10549 | 215 if (isnumeric (m)) # is M a matrix? |
216 if (isempty (m)) # if M is empty, use no precond | |
5837 | 217 p = r; |
10549 | 218 else # otherwise, apply the precond |
8507 | 219 p = m \ r; |
5837 | 220 endif |
10549 | 221 else # then M should be a function! |
8507 | 222 p = feval (m, r, varargin{:}); |
5837 | 223 endif |
224 | |
225 iter = 2; | |
226 | |
227 b_bot_old = 1; | |
5838 | 228 q_old = p_old = s_old = zeros (size (x)); |
5837 | 229 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
230 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
|
231 q = A * p; |
10549 | 232 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
|
233 q = feval (A, p, varargin{:}); |
5837 | 234 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
235 |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
236 resvec(1) = abs (norm (r)); |
5837 | 237 |
238 ## iteration | |
5838 | 239 while (resvec(iter-1) > tol*resvec(1) && iter < maxit) |
240 | |
10549 | 241 if (isnumeric (m)) # is M a matrix? |
242 if (isempty (m)) # if M is empty, use no precond | |
243 s = q; | |
244 else # otherwise, apply the precond | |
245 s = m \ q; | |
5837 | 246 endif |
10549 | 247 else # then M should be a function! |
8507 | 248 s = feval (m, q, varargin{:}); |
5837 | 249 endif |
5838 | 250 b_top = r' * s; |
251 b_bot = q' * s; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
252 |
5837 | 253 if (b_bot == 0.0) |
254 breakdown = true; | |
255 break; | |
256 endif | |
5838 | 257 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
|
258 |
5838 | 259 x += lambda*p; |
260 r -= lambda*q; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
261 |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
262 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
|
263 t = A*s; |
10549 | 264 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
|
265 t = feval (A, s, varargin{:}); |
5837 | 266 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11563
diff
changeset
|
267 |
5838 | 268 alpha0 = (t'*s) / b_bot; |
269 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
|
270 |
5838 | 271 p_temp = p; |
272 q_temp = q; | |
273 | |
5837 | 274 p = s - alpha0*p - alpha1*p_old; |
275 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
|
276 |
5837 | 277 s_old = s; |
278 p_old = p_temp; | |
279 q_old = q_temp; | |
280 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
|
281 |
5838 | 282 resvec(iter) = abs (norm (r)); |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20727
diff
changeset
|
283 iter += 1; |
5837 | 284 endwhile |
285 | |
286 flag = 0; | |
5838 | 287 relres = resvec(iter-1) ./ resvec(1); |
288 iter -= 2; | |
289 if (iter >= maxit-2) | |
5837 | 290 flag = 1; |
291 if (nargout < 2) | |
6498 | 292 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
|
293 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
|
294 1.0/relres); |
5837 | 295 endif |
5838 | 296 elseif (nargout < 2 && ! breakdown) |
6498 | 297 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
|
298 fprintf (stderr, "pcr: the initial residual norm was reduced %g times\n", |
10549 | 299 1.0 / relres); |
5837 | 300 endif |
5838 | 301 |
5837 | 302 if (breakdown) |
303 flag = 3; | |
304 if (nargout < 2) | |
7001 | 305 warning ("pcr: breakdown occurred:\n"); |
6498 | 306 warning ("system matrix singular or preconditioner indefinite?\n"); |
5837 | 307 endif |
308 endif | |
309 | |
310 endfunction | |
311 | |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
312 |
5837 | 313 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
314 %! ## Simplest usage of PCR (see also 'help pcr') |
5837 | 315 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
316 %! N = 20; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
317 %! 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
|
318 %! 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
|
319 %! x = pcr (A,b); |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
320 %! printf ("The solution relative error is %g\n", norm (x-y) / norm (y)); |
5837 | 321 %! |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
322 %! ## 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
|
323 %! ## 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
|
324 %! ## 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
|
325 |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
326 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
327 %! ## Full output from PCR |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
328 %! ## We use this output to plot the convergence history |
5837 | 329 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
330 %! N = 20; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
331 %! 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
|
332 %! 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
|
333 %! [x, flag, relres, iter, resvec] = pcr (A,b); |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
334 %! 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
|
335 %! clf; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
336 %! title ("Convergence history"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
337 %! 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
|
338 %! 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
|
339 |
5837 | 340 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
341 %! ## Full output from PCR |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
342 %! ## 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
|
343 %! ## strongly negative eigenvalue |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
344 %! ## 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
|
345 %! ## and that's why PCR WILL have problems |
5837 | 346 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
347 %! N = 10; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
348 %! 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
|
349 %! X = A \ b; # X is the true solution |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
350 %! 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
|
351 %! [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
|
352 %! if (flag == 3) |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
353 %! 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
|
354 %! endif |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
355 %! clf; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
356 %! title ("Convergence history"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
357 %! 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
|
358 %! 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
|
359 |
5837 | 360 %!demo |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
361 %! ## Full output from PCR |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
362 %! ## 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
|
363 %! ## 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
|
364 %! ## 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
|
365 %! ## 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
|
366 %! ## which is the diagonal of A. |
5837 | 367 %! |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
368 %! ## Note that we use here indefinite preconditioners! |
5837 | 369 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
370 %! N = 100; |
25755
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
371 %! ## Form 1-D Laplacian matrix |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
372 %! A = 2 * eye (N,N); |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
373 %! A(2:(N+1):end) = -1; |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
374 %! 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
|
375 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
376 %! 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
|
377 %! 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
|
378 %! 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
|
379 %! maxit = 80; |
21634
96518f623c91
Backed out changeset dcf8922b724b
Mike Miller <mtmiller@octave.org>
parents:
21633
diff
changeset
|
380 %! 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
|
381 %! ## No preconditioner: the convergence is very slow! |
5837 | 382 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
383 %! [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
|
384 %! clf; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
385 %! title ("Convergence history"); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
386 %! 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
|
387 %! semilogy ([0:iter], resvec, "o-g;NO preconditioning: absolute residual;"); |
5837 | 388 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
389 %! pause (1); |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
390 %! ## Test Jacobi preconditioner: it will not help much!!! |
5837 | 391 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
392 %! 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
|
393 %! [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
|
394 %! hold on; |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
395 %! semilogy ([0:iter],resvec,"o-r;JACOBI preconditioner: absolute residual;"); |
5837 | 396 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
397 %! pause (1); |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
398 %! ## 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
|
399 %! ## 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
|
400 %! |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
401 %! 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
|
402 %! 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
|
403 %! 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
|
404 %! endfor |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
405 %! 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
|
406 %! [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
|
407 %! 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
|
408 %! 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
|
409 |
5837 | 410 %!test |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
411 %! ## solve small indefinite diagonal system |
5837 | 412 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
413 %! 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
|
414 %! 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
|
415 %! 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
|
416 %! [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
|
417 %! 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
|
418 %! 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
|
419 |
5837 | 420 %!test |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
421 %! ## 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
|
422 %! ## should perform max allowable default number of iterations |
5837 | 423 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
424 %! N = 100; |
25755
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
425 %! ## Form 1-D Laplacian matrix |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
426 %! A = 2 * eye (N,N); |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
427 %! A(2:(N+1):end) = -1; |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
428 %! 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
|
429 %! 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
|
430 %! 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
|
431 %! [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
|
432 %! 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
|
433 %! 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
|
434 %! 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
|
435 |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
436 %!test |
17336
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
437 %! ## solve tridiagonal system with "perfect" preconditioner |
b81b9d079515
Use '##' for comments which stand alone on a line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
438 %! ## converges in one iteration |
5837 | 439 %! |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
440 %! N = 100; |
25755
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
441 %! ## Form 1-D Laplacian matrix |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
442 %! A = 2 * eye (N,N); |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
443 %! A(2:(N+1):end) = -1; |
ddf8a03beffb
Calculate Laplacian matrix in BIST tests without for loop.
Rik <rik@octave.org>
parents:
25472
diff
changeset
|
444 %! 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
|
445 %! 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
|
446 %! 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
|
447 %! [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
|
448 %! 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
|
449 %! 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
|
450 %! 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
|
451 %! assert (iter, 1); # should converge in one iteration |