File: gmres.m

package info (click to toggle)
octave 3.8.2-4
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 84,396 kB
  • ctags: 45,547
  • sloc: cpp: 293,356; ansic: 42,041; fortran: 23,669; sh: 13,629; objc: 7,890; yacc: 7,093; lex: 3,442; java: 2,125; makefile: 1,589; perl: 1,009; awk: 974; xml: 34
file content (237 lines) | stat: -rw-r--r-- 7,093 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
## Copyright (C) 2009-2013 Carlo de Falco
##
## This file is part of Octave.
##
## Octave is free software; you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by the
## Free Software Foundation; either version 3 of the License, or (at your
## option) any later version.
##
## Octave is distributed in the hope that it will be useful, but WITHOUT
## ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
## FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
## for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn  {Function File} {@var{x} =} gmres (@var{A}, @var{b}, @var{m}, @var{rtol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0})
## @deftypefnx {Function File} {@var{x} =} gmres (@var{A}, @var{b}, @var{m}, @var{rtol}, @var{maxit}, @var{P})
## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} gmres (@dots{})
## Solve @code{A x = b} using the Preconditioned GMRES iterative method
## with restart, a.k.a. PGMRES(m).
##
## @itemize @minus
## @item @var{rtol} is the relative tolerance,
## if not given or set to [] the default value 1e-6 is used.
##
## @item @var{maxit} is the maximum number of outer iterations,
## if not given or set to [] the default value
## @code{min (10, numel (b) / restart)} is used.
##
## @item @var{x0} is the initial guess,
## if not given or set to [] the default value @code{zeros (size (b))} is used.
##
## @item @var{m} is the restart parameter,
## if not given or set to [] the default value @code{numel (b)} is used.
## @end itemize
##
## Argument @var{A} can be passed as a matrix, function handle, or
## inline function @code{f} such that @code{f(x) = A*x}.
##
## The preconditioner @var{P} is given as @code{P = M1 * M2}.
## Both @var{M1} and @var{M2} can be passed as a matrix, function handle, or
## inline function @code{g} such that @code{g(x) = M1\x} or @code{g(x) = M2\x}.
##
## Besides the vector @var{x}, additional outputs are:
##
## @itemize @minus
## @item @var{flag} indicates the exit status:
##
## @table @asis
## @item 0 : iteration converged to within the specified tolerance
##
## @item 1 : maximum number of iterations exceeded
##
## @item 2 : unused, but skipped for compatibility
##
## @item 3 : algorithm reached stagnation (no change between iterations)
## @end table
##
## @item @var{relres} is the final value of the relative residual.
##
## @item @var{iter} is a vector containing the number of outer iterations and
## total iterations performed.
##
## @item @var{resvec} is a vector containing the relative residual at each
## iteration.
## @end itemize
##
## @seealso{bicg, bicgstab, cgs, pcg}
## @end deftypefn

function [x, flag, relres, it, resvec] = gmres (A, b, restart, rtol, maxit, M1, M2, x0)

  if (nargin < 2 || nargin > 8)
    print_usage ();
  endif

  if (ischar (A))
    Ax = str2func (A);
  elseif (ismatrix (A))
    Ax = @(x) A*x;
  elseif (isa (A, "function_handle"))
    Ax = A;
  else
    error ("gmres: A must be a function or matrix");
  endif

  if (nargin < 3 || isempty (restart))
    restart = rows (b);
  endif

  if (nargin < 4 || isempty (rtol))
    rtol = 1e-6;
  endif

  if (nargin < 5 || isempty (maxit))
    maxit = min (rows (b)/restart, 10);
  endif

  if (nargin < 6 || isempty (M1))
    M1m1x = @(x) x;
  elseif (ischar (M1))
    M1m1x = str2func (M1);
  elseif (ismatrix (M1))
    M1m1x = @(x) M1 \ x;
  elseif (isa (M1, "function_handle"))
    M1m1x = M1;
  else
    error ("gmres: preconditioner M1 must be a function or matrix");
  endif

  if (nargin < 7 || isempty (M2))
    M2m1x = @(x) x;
  elseif (ischar (M2))
    M2m1x = str2func (M2);
  elseif (ismatrix (M2))
    M2m1x = @(x) M2 \ x;
  elseif (isa (M2, "function_handle"))
    M2m1x = M2;
  else
    error ("gmres: preconditioner M2 must be a function or matrix");
  endif

  Pm1x = @(x) M2m1x (M1m1x (x));

  if (nargin < 8 || isempty (x0))
    x0 = zeros (size (b));
  endif

  x_old = x0;
  x = x_old;
  prec_res = Pm1x (b - Ax (x_old));
  presn = norm (prec_res, 2);

  B = zeros (restart + 1, 1);
  V = zeros (rows (x), restart);
  H = zeros (restart + 1, restart);

  ## begin loop
  iter = 1;
  restart_it  = restart + 1;
  resvec      = zeros (maxit, 1);
  resvec(1)   = presn;
  prec_b_norm = norm (Pm1x (b), 2);
  flag        = 1;  # Default flag is maximum # of iterations exceeded

  while (iter <= maxit * restart && presn > rtol * prec_b_norm)

    ## restart
    if (restart_it > restart)
      restart_it = 1;
      x_old = x;
      prec_res = Pm1x (b - Ax (x_old));
      presn = norm (prec_res, 2);
      B(1) = presn;
      H(:) = 0;
      V(:, 1) = prec_res / presn;
    endif

    ## basic iteration
    tmp = Pm1x (Ax (V(:, restart_it)));
    [V(:,restart_it+1), H(1:restart_it+1, restart_it)] = ...
        mgorth (tmp, V(:,1:restart_it));

    Y = (H(1:restart_it+1, 1:restart_it) \ B (1:restart_it+1));

    little_res = B(1:restart_it+1) - ...
        H(1:restart_it+1, 1:restart_it) * Y(1:restart_it);

    presn = norm (little_res, 2);

    x = x_old + V(:, 1:restart_it) * Y(1:restart_it);

    resvec(iter+1) = presn;
    if (norm (x - x_old, inf) <= eps)
      flag = 3;  # Stagnation: no change between iterations
      break;
    endif

    restart_it++ ;
    iter++;
  endwhile

  if (nargout > 1)
    ## Calculate extra outputs as requested
    relres = presn / prec_b_norm;
    if (relres <= rtol)
      flag = 0;  # Converged to solution within tolerance
    endif

    it = [floor(iter/restart), restart_it-1];
  endif

endfunction


%!demo
%! dim = 20;
%! A = spdiags ([-ones(dim,1) 2*ones(dim,1) ones(dim,1)], [-1:1], dim, dim);
%! b = ones (dim, 1);
%! [x, flag, relres, iter, resvec] = gmres (A, b, 10, 1e-10, dim, @(x) x ./ diag (A), [], b)

%!shared A, b, dim
%! dim = 100;
%!test
%! A = spdiags ([-ones(dim,1) 2*ones(dim,1) ones(dim,1)], [-1:1], dim, dim);
%! b = ones (dim, 1);
%! x = gmres (A, b, 10, 1e-10, dim, @(x) x ./ diag (A), [], b);
%! assert (x, A\b, 1e-9*norm (x, Inf));
%!
%!test
%! x = gmres (A, b, dim, 1e-10, 1e4, @(x) diag (diag (A)) \ x, [], b);
%! assert(x, A\b, 1e-7*norm (x, Inf));
%!
%!test
%! 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);
%! A = A'*A;
%! b = rand (dim, 1);
%! [x, resvec] = gmres (@(x) A*x, b, dim, 1e-10, dim, @(x) x./diag (A), [], []);
%! assert (x, A\b, 1e-9*norm (x, Inf));
%! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) diag (diag (A)) \ x, [], []);
%! assert (x, A\b, 1e-9*norm (x, Inf));
%!test
%! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) x ./ diag (A), [], []);
%! assert (x, A\b, 1e-7*norm (x, Inf));


%!error gmres (1)
%!error gmres (1,2,3,4,5,6,7,8,9)
%!error <A must be> gmres ({1},2)
%!error <A must be a function or matrix> gmres ({1},2)
%!error <M1 must be a function or matrix> gmres (1,2,3,4,5,{6})
%!error <M2 must be a function or matrix> gmres (1,2,3,4,5,6,{7})