annotate scripts/linear-algebra/gmres.m @ 12541:dd2c70b30f28

Add tests for ifftshift.m
author Robert T. Short <octave@phaselockedsystems.com.com>
date Sat, 26 Mar 2011 06:50:12 -0700
parents e742720c5e71
children f96b9b9f141b
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
12509
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
1 ## Copyright (C) 2009-2011 Carlo de Falco
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
2 ##
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.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
4 ##
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.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
9 ##
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.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
14 ##
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/>.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
18
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{})
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
23 ## Solve @code{A x = b} using the Preconditioned GMRES iterative method
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
24 ## with restart, a.k.a. PGMRES(m).
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
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
27 ## @item @var{rtol} is the relative tolerance,
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 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
30 ## @item @var{maxit} is the maximum number of outer iterations,
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
31 ## if not given or set to [] the default value
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
32 ## @code{min (10, numel (b) / restart)} is used.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
33 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
34 ## @item @var{x0} is the initial guess,
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
35 ## if not given or set to [] the default value @code{zeros(size (b))} is used.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
36 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
37 ## @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
38 ## 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
39 ## @end itemize
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
40 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
41 ## Argument @var{A} can be passed as a matrix, function handle, or
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
42 ## inline function @code{f} such that @code{f(x) = A*x}.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
43 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
44 ## The preconditioner @var{P} is given as @code{P = M1 * M2}.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
45 ## Both @var{M1} and @var{M2} can be passed as a matrix, function handle, or
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
46 ## inline function @code{g} such that @code{g(x) = M1\x} or @code{g(x) = M2\x}.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
47 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
48 ## 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
49 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
50 ## @itemize @minus
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
51 ## @item @var{flag} indicates the exit status:
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
52 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
53 ## @table @asis
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
54 ## @item 0 : iteration converged to within the specified tolerance
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
55 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
56 ## @item 1 : maximum number of iterations exceeded
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
57 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
58 ## @item 2 : unused, but skipped for compatibility
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
59 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
60 ## @item 3 : algorithm reached stagnation
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
61 ## @end table
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
62 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
63 ## @item @var{relres} is the final value of the relative residual.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
64 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
65 ## @item @var{iter} is a vector containing the number of outer iterations and
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
66 ## total iterations performed.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
67 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
68 ## @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
69 ## iteration.
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
70 ## @end itemize
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
71 ##
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
72 ## @seealso{pcg, cgs, bigcstab}
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
73 ## @end deftypefn
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
74
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
75 function [x, flag, prec_res_norm, itcnt] = gmres (A, b, restart, rtol, maxit, M1, M2, x0)
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
76
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
77 if (nargin < 2 || nargin > 8)
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
78 print_usage ();
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
79 end
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
80
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
81 if (ischar (A))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
82 Ax = str2func (A);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
83 elseif (ismatrix (A))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
84 Ax = @(x) A*x;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
85 elseif (isa (A, "function_handle"))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
86 Ax = A;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
87 else
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
88 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
89 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
90
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
91 if (nargin < 3 || isempty (restart))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
92 restart = rows (b);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
93 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
94
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
95 if (nargin < 4 || isempty (rtol))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
96 rtol = 1e-6;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
97 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
98
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
99 if (nargin < 5 || isempty (maxit))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
100 maxit = min (rows (b)/restart, 10);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
101 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
102
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
103 if (nargin < 6 || isempty (M1))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
104 M1m1x = @(x) x;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
105 elseif (ischar (M1))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
106 M1m1x = str2func (M1);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
107 elseif (ismatrix (M1))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
108 M1m1x = @(x) M1 \ x;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
109 elseif (isa (M1, "function_handle"))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
110 M1m1x = M1;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
111 else
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
112 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
113 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
114
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
115 if (nargin < 7 || isempty (M2))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
116 M2m1x = @(x) x;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
117 elseif (ischar (M2))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
118 M2m1x = str2func (M2);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
119 elseif (ismatrix (M2))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
120 M2m1x = @(x) M2 \ x;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
121 elseif (isa (M2, "function_handle"))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
122 M2m1x = M2;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
123 else
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
124 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
125 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
126
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
127 Pm1x = @(x) M2m1x (M1m1x (x));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
128
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
129 if (nargin < 8 || isempty (x0))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
130 x0 = zeros (size (b));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
131 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
132
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
133 x_old = x0;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
134 x = x_old;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
135 prec_res = Pm1x (b - Ax (x_old));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
136 prec_res_norm = norm (prec_res, 2);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
137
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
138 B = zeros (restart + 1, 1);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
139 V = zeros (rows (x), restart);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
140 H = zeros (restart + 1, restart);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
141
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
142 ## begin loop
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
143 iter = 1;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
144 restart_it = restart + 1;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
145 resids = zeros (maxit, 1);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
146 resids(1) = prec_res_norm;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
147 prec_b_norm = norm (Pm1x (b), 2);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
148 flag = 1;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
149
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
150 while ((iter <= maxit * restart) && (prec_res_norm > rtol * prec_b_norm))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
151
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
152 ## restart
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
153 if (restart_it > restart)
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
154 restart_it = 1;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
155 x_old = x;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
156 prec_res = Pm1x (b - Ax (x_old));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
157 prec_res_norm = norm (prec_res, 2);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
158 B(1) = prec_res_norm;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
159 H(:) = 0;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
160 V(:, 1) = prec_res / prec_res_norm;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
161 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
162
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
163 ## basic iteration
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
164 tmp = Pm1x (Ax (V(:, restart_it)));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
165 [V(:,restart_it+1), H(1:restart_it+1, restart_it)] = mgorth (tmp, V(:,1:restart_it));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
166
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
167 Y = (H(1:restart_it+1, 1:restart_it) \ B (1:restart_it+1));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
168
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
169 little_res = B(1:restart_it+1) - H(1:restart_it+1, 1:restart_it) * Y(1:restart_it);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
170 prec_res_norm = norm (little_res, 2);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
171
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
172 x = x_old + V(:, 1:restart_it) * Y(1:restart_it);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
173
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
174 resids(iter) = prec_res_norm ;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
175 if (norm (x - x_old, inf) <= eps)
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
176 flag = 3;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
177 break
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
178 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
179
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
180 restart_it++ ; iter++;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
181 endwhile
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
182
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
183 if (prec_res_norm > rtol * prec_b_norm)
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
184 flag = 0;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
185 endif
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
186
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
187 resids = resids(1:iter-1);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
188 itcnt = [floor(maxit/restart), rem(maxit, restart)];
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
189 endfunction
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
190
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
191
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
192 %!shared A, b, dim
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
193 %! dim = 100;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
194 %!test
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
195 %! A = spdiags ([-ones(dim,1) 2*ones(dim,1) ones(dim,1)], [-1:1], dim, dim);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
196 %! b = ones(dim, 1);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
197 %! x = gmres (A, b, 10, 1e-10, dim, @(x) x./diag(A), [], b);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
198 %! assert(x, A\b, 1e-9*norm(x,inf));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
199 %!
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
200 %!test
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
201 %! x = gmres (A, b, dim, 1e-10, 1e4, @(x) diag(diag(A))\x, [], b);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
202 %! assert(x, A\b, 1e-7*norm(x,inf));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
203 %!
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
204 %!test
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
205 %! 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
206 %! A = A'*A;
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
207 %! b = rand (dim, 1);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
208 %! [x, resids] = gmres (@(x) A*x, b, dim, 1e-10, dim, @(x) x./diag(A), [], []);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
209 %! assert(x, A\b, 1e-9*norm(x,inf))
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
210 %! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) diag(diag(A))\x, [], []);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
211 %! assert(x, A\b, 1e-9*norm(x,inf));
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
212 %!test
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
213 %! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) x./diag(A), [], []);
e742720c5e71 Add mgorth and gmres functions from Octave Forge
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff changeset
214 %! assert(x, A\b, 1e-7*norm(x,inf));