Mercurial > hg > octave-nkf
annotate scripts/sparse/gmres.m @ 20761:b7ac1e94266e
maint: Further clean up of functions in ode/private dir.
* AbsRel_Norm.m, fuzzy_compare.m, integrate_adaptive.m, integrate_const.m,
integrate_n_steps.m, ode_struct_value_check.m, odepkg_event_handle.m,
odepkg_structure_check.m, runge_kutta_45_dorpri.m:
Place latest copyright first in file.
Use two spaces before beginning single-line comment.
Use parentheses around variable to be tested in switch stmt.
Use space between function name and opening parenthesis.
author | Rik <rik@octave.org> |
---|---|
date | Mon, 05 Oct 2015 12:03:16 -0700 |
parents | df437a52bcaf |
children |
rev | line source |
---|---|
19898
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
1 ## Copyright (C) 2009-2015 Carlo de Falco |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
2 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
3 ## This file is part of Octave. |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
4 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
5 ## Octave is free software; you can redistribute it and/or modify it |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
6 ## under the terms of the GNU General Public License as published by the |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
7 ## Free Software Foundation; either version 3 of the License, or (at your |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
8 ## option) any later version. |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
9 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
10 ## Octave is distributed in the hope that it will be useful, but WITHOUT |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
11 ## ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
12 ## FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
13 ## for more details. |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
14 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
15 ## You should have received a copy of the GNU General Public License |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
16 ## along with Octave; see the file COPYING. If not, see |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
17 ## <http://www.gnu.org/licenses/>. |
12761
13bcd62824a7
doc: Add documentation for gmres, rectangle to manual
Rik <octave@nomad.inbox5.com>
parents:
12642
diff
changeset
|
18 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
19 ## -*- texinfo -*- |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
20 ## @deftypefn {Function File} {@var{x} =} gmres (@var{A}, @var{b}, @var{m}, @var{rtol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0}) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
21 ## @deftypefnx {Function File} {@var{x} =} gmres (@var{A}, @var{b}, @var{m}, @var{rtol}, @var{maxit}, @var{P}) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
22 ## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} gmres (@dots{}) |
20374
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
23 ## Solve @code{A x = b} using the Preconditioned GMRES iterative method with |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
24 ## restart, a.k.a. PGMRES(m). |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
25 ## |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
26 ## @itemize @minus |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
12509
diff
changeset
|
27 ## @item @var{rtol} is the relative tolerance, |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
28 ## if not given or set to [] the default value 1e-6 is used. |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
29 ## |
20374
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
30 ## @item @var{maxit} is the maximum number of outer iterations, if not given or |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
31 ## set to [] the default value @code{min (10, numel (b) / restart)} is used. |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
32 ## |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
33 ## @item @var{x0} is the initial guess, |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
34 ## if not given or set to [] the default value @code{zeros (size (b))} is used. |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
35 ## |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
36 ## @item @var{m} is the restart parameter, |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
37 ## if not given or set to [] the default value @code{numel (b)} is used. |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
38 ## @end itemize |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
39 ## |
20374
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
40 ## Argument @var{A} can be passed as a matrix, function handle, or inline |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
41 ## function @code{f} such that @code{f(x) = A*x}. |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
42 ## |
20374
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
43 ## The preconditioner @var{P} is given as @code{P = M1 * M2}. Both @var{M1} |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
44 ## and @var{M2} can be passed as a matrix, function handle, or inline function |
df437a52bcaf
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
20059
diff
changeset
|
45 ## @code{g} such that @code{g(x) = M1\x} or @code{g(x) = M2\x}. |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
46 ## |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
47 ## Besides the vector @var{x}, additional outputs are: |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
48 ## |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
49 ## @itemize @minus |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
50 ## @item @var{flag} indicates the exit status: |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
51 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
52 ## @table @asis |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
53 ## @item 0 : iteration converged to within the specified tolerance |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13861
diff
changeset
|
54 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
55 ## @item 1 : maximum number of iterations exceeded |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13861
diff
changeset
|
56 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
57 ## @item 2 : unused, but skipped for compatibility |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13861
diff
changeset
|
58 ## |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
59 ## @item 3 : algorithm reached stagnation (no change between iterations) |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
60 ## @end table |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13861
diff
changeset
|
61 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
62 ## @item @var{relres} is the final value of the relative residual. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13861
diff
changeset
|
63 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
64 ## @item @var{iter} is a vector containing the number of outer iterations and |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
12509
diff
changeset
|
65 ## total iterations performed. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13861
diff
changeset
|
66 ## |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
67 ## @item @var{resvec} is a vector containing the relative residual at each |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
68 ## iteration. |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
69 ## @end itemize |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
70 ## |
20059 | 71 ## @seealso{bicg, bicgstab, cgs, pcg, pcr, qmr} |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
72 ## @end deftypefn |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
73 |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
74 function [x, flag, relres, it, resvec] = gmres (A, b, restart, rtol, maxit, M1, M2, x0) |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
75 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
76 if (nargin < 2 || nargin > 8) |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
77 print_usage (); |
13174
bd2cd4fd3edf
maint: use specific endif, endfor tokens instead of simple end
John W. Eaton <jwe@octave.org>
parents:
12761
diff
changeset
|
78 endif |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
79 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
80 if (ischar (A)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
81 Ax = str2func (A); |
19901
00e31f316a3a
Fix Matlab incompatibility of "ismatrix" (bug #42422).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
19898
diff
changeset
|
82 elseif (isnumeric (A) && ismatrix (A)) |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
83 Ax = @(x) A*x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
84 elseif (isa (A, "function_handle")) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
85 Ax = A; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
86 else |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
87 error ("gmres: A must be a function or matrix"); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
88 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
89 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
90 if (nargin < 3 || isempty (restart)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
91 restart = rows (b); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
92 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
93 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
94 if (nargin < 4 || isempty (rtol)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
95 rtol = 1e-6; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
96 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
97 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
98 if (nargin < 5 || isempty (maxit)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
99 maxit = min (rows (b)/restart, 10); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
100 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
101 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
102 if (nargin < 6 || isempty (M1)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
103 M1m1x = @(x) x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
104 elseif (ischar (M1)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
105 M1m1x = str2func (M1); |
19901
00e31f316a3a
Fix Matlab incompatibility of "ismatrix" (bug #42422).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
19898
diff
changeset
|
106 elseif (isnumeric (M1) && ismatrix (M1)) |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
107 M1m1x = @(x) M1 \ x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
108 elseif (isa (M1, "function_handle")) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
109 M1m1x = M1; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
110 else |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
111 error ("gmres: preconditioner M1 must be a function or matrix"); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
112 endif |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
113 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
114 if (nargin < 7 || isempty (M2)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
115 M2m1x = @(x) x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
116 elseif (ischar (M2)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
117 M2m1x = str2func (M2); |
19901
00e31f316a3a
Fix Matlab incompatibility of "ismatrix" (bug #42422).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
19898
diff
changeset
|
118 elseif (isnumeric (M2) && ismatrix (M2)) |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
119 M2m1x = @(x) M2 \ x; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
120 elseif (isa (M2, "function_handle")) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
121 M2m1x = M2; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
122 else |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
123 error ("gmres: preconditioner M2 must be a function or matrix"); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
124 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
125 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
126 Pm1x = @(x) M2m1x (M1m1x (x)); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
127 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
128 if (nargin < 8 || isempty (x0)) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
129 x0 = zeros (size (b)); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
130 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
131 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
132 x_old = x0; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
133 x = x_old; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
134 prec_res = Pm1x (b - Ax (x_old)); |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
135 presn = norm (prec_res, 2); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
136 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
137 B = zeros (restart + 1, 1); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
138 V = zeros (rows (x), restart); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
139 H = zeros (restart + 1, restart); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
140 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
141 ## begin loop |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
142 iter = 1; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
143 restart_it = restart + 1; |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
144 resvec = zeros (maxit, 1); |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
145 resvec(1) = presn; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
146 prec_b_norm = norm (Pm1x (b), 2); |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
147 flag = 1; # Default flag is maximum # of iterations exceeded |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
148 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
149 while (iter <= maxit * restart && presn > rtol * prec_b_norm) |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
150 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
151 ## restart |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
152 if (restart_it > restart) |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
153 restart_it = 1; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
154 x_old = x; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
155 prec_res = Pm1x (b - Ax (x_old)); |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
156 presn = norm (prec_res, 2); |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
157 B(1) = presn; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
158 H(:) = 0; |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
159 V(:, 1) = prec_res / presn; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
160 endif |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
161 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
162 ## basic iteration |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
163 tmp = Pm1x (Ax (V(:, restart_it))); |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
164 [V(:,restart_it+1), H(1:restart_it+1, restart_it)] = ... |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
165 mgorth (tmp, V(:,1:restart_it)); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
166 |
20038
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19901
diff
changeset
|
167 Y = (H(1:restart_it+1, 1:restart_it) \ B(1:restart_it+1)); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
168 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
169 little_res = B(1:restart_it+1) - ... |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
170 H(1:restart_it+1, 1:restart_it) * Y(1:restart_it); |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
171 |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
172 presn = norm (little_res, 2); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
173 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
174 x = x_old + V(:, 1:restart_it) * Y(1:restart_it); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
175 |
14822
93e5ade3fda4
Improve matlab compatibility of gmres
Carlo de Falco <cdf@users.sourceforge.net>
parents:
14773
diff
changeset
|
176 resvec(iter+1) = presn; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
177 if (norm (x - x_old, inf) <= eps) |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
178 flag = 3; # Stagnation: no change between iterations |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
179 break; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
180 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
181 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
182 restart_it++ ; |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
12763
diff
changeset
|
183 iter++; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
184 endwhile |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
185 |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
186 if (nargout > 1) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
187 ## Calculate extra outputs as requested |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
188 relres = presn / prec_b_norm; |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
189 if (relres <= rtol) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
190 flag = 0; # Converged to solution within tolerance |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
191 endif |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
192 |
14822
93e5ade3fda4
Improve matlab compatibility of gmres
Carlo de Falco <cdf@users.sourceforge.net>
parents:
14773
diff
changeset
|
193 it = [floor(iter/restart), restart_it-1]; |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
194 endif |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
195 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
196 endfunction |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
197 |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
198 |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
199 %!demo |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
200 %! dim = 20; |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
201 %! A = spdiags ([-ones(dim,1) 2*ones(dim,1) ones(dim,1)], [-1:1], dim, dim); |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
202 %! b = ones (dim, 1); |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
203 %! [x, flag, relres, iter, resvec] = gmres (A, b, 10, 1e-10, dim, @(x) x ./ diag (A), [], b) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
204 |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
205 %!shared A, b, dim |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
206 %! dim = 100; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
207 %!test |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
208 %! A = spdiags ([-ones(dim,1) 2*ones(dim,1) ones(dim,1)], [-1:1], dim, dim); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
209 %! b = ones (dim, 1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
210 %! x = gmres (A, b, 10, 1e-10, dim, @(x) x ./ diag (A), [], b); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
211 %! assert (x, A\b, 1e-9*norm (x, Inf)); |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
212 %! |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
213 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
214 %! x = gmres (A, b, dim, 1e-10, 1e4, @(x) diag (diag (A)) \ x, [], b); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
215 %! assert(x, A\b, 1e-7*norm (x, Inf)); |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
216 %! |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
217 %!test |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
218 %! A = spdiags ([[1./(2:2:2*(dim-1)) 0]; 1./(1:2:2*dim-1); [0 1./(2:2:2*(dim-1))]]', -1:1, dim, dim); |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
219 %! A = A'*A; |
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
220 %! b = rand (dim, 1); |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
221 %! [x, resvec] = gmres (@(x) A*x, b, dim, 1e-10, dim, @(x) x./diag (A), [], []); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
222 %! assert (x, A\b, 1e-9*norm (x, Inf)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
223 %! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) diag (diag (A)) \ x, [], []); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
224 %! assert (x, A\b, 1e-9*norm (x, Inf)); |
12509
e742720c5e71
Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
225 %!test |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14822
diff
changeset
|
226 %! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) x ./ diag (A), [], []); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
227 %! assert (x, A\b, 1e-7*norm (x, Inf)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
228 |
14773
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
229 |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
230 %!error gmres (1) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
231 %!error gmres (1,2,3,4,5,6,7,8,9) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
232 %!error <A must be> gmres ({1},2) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
233 %!error <A must be a function or matrix> gmres ({1},2) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
234 %!error <M1 must be a function or matrix> gmres (1,2,3,4,5,{6}) |
5c269b73f467
gmres.m: Overhaul function. Fix bug #36568.
Rik <octave@nomad.inbox5.com>
parents:
14366
diff
changeset
|
235 %!error <M2 must be a function or matrix> gmres (1,2,3,4,5,6,{7}) |
17338
1c89599167a6
maint: End m-files with 1 blank line.
Rik <rik@octave.org>
parents:
14868
diff
changeset
|
236 |