annotate scripts/sparse/pcr.m @ 19697:4197fc428c7d

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