File: qmr.m

package info (click to toggle)
octave 9.4.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 144,300 kB
  • sloc: cpp: 332,784; ansic: 77,239; fortran: 20,963; objc: 9,396; sh: 8,213; yacc: 4,925; lex: 4,389; perl: 1,544; java: 1,366; awk: 1,259; makefile: 648; xml: 189
file content (328 lines) | stat: -rw-r--r-- 9,729 bytes parent folder | download
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
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
########################################################################
##
## Copyright (C) 2014-2024 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## 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
## <https://www.gnu.org/licenses/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn  {} {@var{x} =} qmr (@var{A}, @var{b}, @var{rtol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0})
## @deftypefnx {} {@var{x} =} qmr (@var{A}, @var{b}, @var{rtol}, @var{maxit}, @var{P})
## @deftypefnx {} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} qmr (@var{A}, @var{b}, @dots{})
## Solve @code{A x = b} using the Quasi-Minimal Residual iterative method
## (without look-ahead).
##
## @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} the maximum number of outer iterations, if not given or
## set to [] the default value @code{min (20, numel (b))} is used.
##
## @item @var{x0} the initial guess, if not given or set to [] the default
## value @code{zeros (size (b))} is used.
## @end itemize
##
## @var{A} can be passed as a matrix or as a function handle or inline
## function @code{f} such that @code{f(x, "notransp") = A*x} and
## @code{f(x, "transp") = 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 or as a function handle or inline
## function @code{g} such that @code{g(x, "notransp") = M1 \ x} or
## @code{g(x, "notransp") = M2 \ x} and @code{g(x, "transp") = M1' \ x} or
## @code{g(x, "transp") = M2' \ x}.
##
## If called with more than one output parameter
##
## @itemize @minus
## @item @var{flag} indicates the exit status:
##
## @itemize @minus
## @item 0: iteration converged to the within the chosen tolerance
##
## @item 1: the maximum number of iterations was reached before convergence
##
## @item 3: the algorithm reached stagnation
## @end itemize
##
## (the value 2 is unused but skipped for compatibility).
##
## @item @var{relres} is the final value of the relative residual.
##
## @item @var{iter} is the number of iterations performed.
##
## @item @var{resvec} is a vector containing the residual norms at each
##       iteration.
## @end itemize
##
## References:
##
## @enumerate
## @item
## @nospell{R. Freund and N. Nachtigal}, @cite{QMR: a quasi-minimal residual
## method for non-Hermitian linear systems}, @nospell{Numerische Mathematik},
## 1991, 60, pp.@: 315--339.
##
## @item
## @nospell{ R. Barrett, M. Berry, T. Chan, J. Demmel, J. Donato, J. Dongarra},
## @nospell{ V. Eijkhour, R. Pozo, C. Romine, and H. van der Vorst},
## @cite{Templates for the solution of linear systems: Building blocks
## for iterative methods}, SIAM, 2nd ed., 1994.
## @end enumerate
##
## @seealso{bicg, bicgstab, cgs, gmres, pcg}
## @end deftypefn

function [x, flag, relres, iter, resvec] = qmr (A, b, rtol, maxit, M1, M2, x0)

  if (nargin >= 2 && isvector (full (b)))

    if (ischar (A))
      fcn = str2func (A);
      Ax  = @(x) feval (fcn, x, "notransp");
      Atx = @(x) feval (fcn, x, "transp");
    elseif (is_function_handle (A))
      Ax  = @(x) feval (A, x, "notransp");
      Atx = @(x) feval (A, x, "transp");
    elseif (isnumeric (A) && issquare (A))
      Ax  = @(x) A  * x;
      Atx = @(x) A' * x;
    else
      error ("qmr: A must be a square matrix or function");
    endif

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

    if (nargin < 4 || isempty (maxit))
      maxit = min (rows (b), 20);
    else
      maxit = fix (maxit);
    endif

    if (nargin < 5 || isempty (M1))
      M1m1x = @(x, ignore) x;
      M1tm1x = M1m1x;
    elseif (ischar (M1))
      fcn = str2func (M1);
      M1m1x  = @(x) feval (fcn, x, "notransp");
      M1tm1x = @(x) feval (fcn, x, "transp");
    elseif (is_function_handle (M1))
      M1m1x  = @(x) feval (M1, x, "notransp");
      M1tm1x = @(x) feval (M1, x, "transp");
    elseif (isnumeric (M1) && ismatrix (M1))
      M1m1x  = @(x) M1  \ x;
      M1tm1x = @(x) M1' \ x;
    else
      error ("qmr: preconditioner M1 must be a function or matrix");
    endif

    if (nargin < 6 || isempty (M2))
      M2m1x = @(x, ignore) x;
      M2tm1x = M2m1x;
    elseif (ischar (M2))
      fcn = str2func (M2);
      M2m1x  = @(x) feval (fcn, x, "notransp");
      M2tm1x = @(x) feval (fcn, x, "transp");
    elseif (is_function_handle (M2))
      M2m1x  = @(x) feval (M2, x, "notransp");
      M2tm1x = @(x) feval (M2, x, "transp");
    elseif (isnumeric (M2) && ismatrix (M2))
      M2m1x  = @(x) M2  \ x;
      M2tm1x = @(x) M2' \ x;
    else
      error ("qmr: preconditioner M2 must be a function or matrix");
    endif

    if (nargin < 7 || isempty (x0))
      x = zeros (size (b));
    else
      x = x0;
    endif

    r = b - Ax (x);

    bnorm = norm (b);
    res0 = norm (r);
    if (nargout > 4)
      resvec(1) = res0;
    endif
    vt = r;

    y = M1m1x (vt);

    rho0 = norm (y);
    wt = r;

    z = M2tm1x (wt);

    xi1 = norm (z);
    gamma0 = 1;
    eta0 = -1;
    flag = 1;
    for iter=1:1:maxit
      ## If rho0 == 0 or xi1 == 0, method fails.
      v = vt / rho0;
      y /= rho0;
      w = wt / xi1;
      z /= xi1;

      delta1 = z' * y;   # If delta1 == 0, method fails.

      yt = M2m1x (y);
      zt = M1tm1x (z);

      if (iter == 1)
        p = yt;
        q = zt;
      else
        p = yt - (xi1*delta1/eps0) * p;
        q = zt - (rho0*delta1/eps0) * q;
      endif
      pt = Ax (p);

      eps0 = q' * pt;          # If eps0 == 0, method fails.
      beta1 = eps0 / delta1;   # If beta1 == 0, method fails.
      vt = pt - beta1 * v;

      y = M1m1x (vt);
      rho1 = norm (y);
      wt = Atx (q) - beta1 * w;
      z = M2tm1x (wt);

      xi1 = norm (z);
      theta1 = rho1 / (gamma0 * abs (beta1));
      gamma1 = 1 / sqrt (1 + theta1^2);   # If gamma1 == 0, method fails.
      eta1 = -eta0 * rho0 * gamma1^2 / (beta1 * gamma0^2);

      if (iter == 1)
        d = eta1 * p;
        s = eta1 * pt;
      else
        d = eta1 * p + (theta0*gamma1)^2 * d;
        s = eta1 * pt + (theta0 * gamma1)^2 * s;
      endif
      x += d;
      r -= s;

      res1 = norm (r) / bnorm;
      if (nargout > 4)
        resvec(iter + 1, 1) = norm (r);
      endif

      if (res1 < rtol)
        ## Convergence achieved.
        flag = 0;
        break;
      elseif (res0 <= res1)
        ## Stagnation encountered.
        flag = 3;
        break;
      endif
      theta0 = theta1;
      eta0 = eta1;
      gamma0 = gamma1;
      rho0 = rho1;
    endfor

    relres = res1;
    if (flag == 1)
      if (nargout < 2)
        printf ("qmr stopped at iteration %i ", iter);
        printf ("without converging to the desired tolerance %e\n", rtol);
        printf ("because the maximum number of iterations was reached. ");
        printf ("The iterate returned (number %i) has ", maxit);
        printf ("relative residual %e\n", res1);
      endif
    elseif (flag == 3)
      if (nargout < 2)
        printf ("qmr stopped at iteration %i ", iter);
        printf (" without converging to the desired tolerance %e\n", rtol);
        printf ("because the method stagnated.\n");
        printf ("The iterate returned (number %i) ", iter);
        printf ("has relative residual %e\n", res1);
      endif
    elseif (nargout < 2)
      printf ("qmr converged at iteration %i ", iter);
      printf ("to a solution with relative residual %e\n", res1);
    endif
  else
    print_usage ();
  endif

endfunction


%!demo
%! % Solve system of A*x=b
%! A = [5 -1 3;-1 2 -2;3 -2 3];
%! b = [7;-1;4];
%! [x, flag, relres, iter, resvec] = qmr (A, b)

%!test
%! n = 100;
%! A = spdiags ([-2*ones(n,1) 4*ones(n,1) -ones(n,1)], -1:1, n, n);
%! b = sum (A, 2);
%! rtol = 1e-8;
%! maxit = 15;
%! M1 = spdiags ([ones(n,1)/(-2) ones(n,1)],-1:0, n, n);
%! M2 = spdiags ([4*ones(n,1) -ones(n,1)], 0:1, n, n);
%! [x, flag, relres, iter, resvec] = qmr (A, b, rtol, maxit, M1, M2);
%! assert (x, ones (size (b)), 1e-7);

%!function y = afcn (x, t, a)
%!  switch (t)
%!    case "notransp"
%!      y = a * x;
%!    case "transp"
%!      y = a' * x;
%!  endswitch
%!endfunction
%!
%!test
%! n = 100;
%! A = spdiags ([-2*ones(n,1) 4*ones(n,1) -ones(n,1)], -1:1, n, n);
%! b = sum (A, 2);
%! rtol = 1e-8;
%! maxit = 15;
%! M1 = spdiags ([ones(n,1)/(-2) ones(n,1)],-1:0, n, n);
%! M2 = spdiags ([4*ones(n,1) -ones(n,1)], 0:1, n, n);
%!
%! [x, flag, relres, iter, resvec] = qmr (@(x, t) afcn (x, t, A),
%!                                         b, rtol, maxit, M1, M2);
%! assert (x, ones (size (b)), 1e-7);

%!test
%! n = 100;
%! rtol = 1e-8;
%! a = sprand (n, n, .1);
%! A = a' * a + 100 * eye (n);
%! b = sum (A, 2);
%! [x, flag, relres, iter, resvec] = qmr (A, b, rtol, [], diag (diag (A)));
%! assert (x, ones (size (b)), 1e-7);

%!test
%! A = [1 + 1i, 1 + 1i; 2 - 1i, 2 + 1i];
%! b = A * [1; 1];
%! [x, flag, relres, iter, resvec] = qmr (A, b);
%! assert (x, [1; 1], 1e-6);