1 ## Copyright (C) 2008-2012 Radek Salac
2 ## Copyright (C) 2012 Carlo de Falco
4 ## This file is part of Octave.
6 ## Octave is free software; you can redistribute it and/or modify it
7 ## under the terms of the GNU General Public License as published by
8 ## the Free Software Foundation; either version 3 of the License, or (at
9 ## your option) any later version.
11 ## Octave is distributed in the hope that it will be useful, but
12 ## WITHOUT ANY WARRANTY; without even the implied warranty of
13 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 ## General Public License for more details.
16 ## You should have received a copy of the GNU General Public License
17 ## along with Octave; see the file COPYING. If not, see
18 ## <http://www.gnu.org/licenses/>.
22 ## @deftypefn {Function File} {@var{x} =} bicgstab (@var{A}, @var{b}, @var{rtol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0})
23 ## @deftypefnx {Function File} {@var{x} =} bicgstab (@var{A}, @var{b}, @var{rtol}, @var{maxit}, @var{P})
24 ## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} bicgstab (@var{A}, @var{b}, @dots{})
25 ## Solve @code{A x = b} using the stabilizied Bi-conjugate gradient iterative
29 ## @item @var{rtol} is the relative tolerance, if not given or set to
30 ## [] the default value 1e-6 is used.
32 ## @item @var{maxit} the maximum number of outer iterations, if not
33 ## given or set to [] the default value @code{min (20, numel (b))} is
36 ## @item @var{x0} the initial guess, if not given or set to [] the
37 ## default value @code{zeros (size (b))} is used.
40 ## @var{A} can be passed as a matrix or as a function handle or
41 ## inline function @code{f} such that @code{f(x) = A*x}.
43 ## The preconditioner @var{P} is given as @code{P = M1 * M2}.
44 ## Both @var{M1} and @var{M2} can be passed as a matrix or as a function
45 ## handle or inline function @code{g} such that @code{g(x) = M1 \ x} or
46 ## @code{g(x) = M2 \ x}.
48 ## If called with more than one output parameter
51 ## @item @var{flag} indicates the exit status:
53 ## @item 0: iteration converged to the within the chosen tolerance
55 ## @item 1: the maximum number of iterations was reached before convergence
57 ## @item 3: the algorithm reached stagnation
59 ## (the value 2 is unused but skipped for compatibility).
61 ## @item @var{relres} is the final value of the relative residual.
63 ## @item @var{iter} is the number of iterations performed.
65 ## @item @var{resvec} is a vector containing the relative residual at each iteration.
68 ## @seealso{bicg, cgs, gmres, pcg}
72 function [x, flag, relres, iter, resvec] = bicgstab (A, b, tol, maxit,
75 if (nargin >= 2 && nargin <= 7 && isvector (full (b)))
81 elseif (isa (A, "function_handle"))
82 Ax = @(x) feval (A, x);
84 error (["bicgstab: first argument is expected " ...
85 "to be a function or a square matrix"]);
88 if (nargin < 3 || isempty (tol))
92 if (nargin < 4 || isempty (maxit))
93 maxit = min (rows (b), 20);
96 if (nargin < 5 || isempty (M1))
99 M1m1x = str2func (M1);
100 elseif (ismatrix (M1))
102 elseif (isa (M1, "function_handle"))
103 M1m1x = @(x) feval (M1, x);
105 error (["bicgstab: preconditioner is " ...
106 "expected to be a function or matrix"]);
109 if (nargin < 6 || isempty (M2))
112 M2m1x = str2func (M2);
113 elseif (ismatrix (M2))
115 elseif (isa (M2, "function_handle"))
116 M2m1x = @(x) feval (M2, x);
118 error (["bicgstab: preconditioner is "...
119 "expected to be a function or matrix"]);
122 precon = @(x) M2m1x (M1m1x (x));
124 if (nargin < 7 || isempty (x0))
125 x0 = zeros (size (b));
128 ## specifies initial estimate x0
130 x = zeros (rows (b), 1);
140 ## Vector of the residual norms for each iteration.
141 resvec = norm(res) / norm_b;
143 ## Default behaviour we don't reach tolerance tol within maxit iterations.
152 beta = (rho_1 / rho_2) * (alpha / omega);
153 p = res + beta * (p - omega * v);
159 alpha = rho_1 / (rr' * v);
165 omega = (t' * s) / (t' * t);
166 x = x + alpha * phat + omega * shat;
170 relres = norm (res) / norm_b;
171 resvec = [resvec; relres];
174 ## We reach tolerance tol within maxit iterations.
177 elseif (resvec(end) == resvec(end - 1))
178 ## The method stagnates.
186 printf ("bicgstab converged at iteration %i ", iter);
187 printf ("to a solution with relative residual %e\n", relres);
189 printf ("bicgstab stopped at iteration %i ", iter);
190 printf ("without converging to the desired tolerance %e\n", tol);
191 printf ("because the method stagnated.\n");
192 printf ("The iterate returned (number %i) ", iter);
193 printf ("has relative residual %e\n", relres);
195 printf ("bicgstab stopped at iteration %i ", iter);
196 printf ("without converging to the desired toleranc %e\n", tol);
197 printf ("because the maximum number of iterations was reached.\n");
198 printf ("The iterate returned (number %i) ", iter);
199 printf ("has relative residual %e\n", relres);
210 %! % Solve system of A*x=b
211 %! A = [5 -1 3;-1 2 -2;3 -2 3]
213 %! [x, flag, relres, iter, resvec] = bicgstab(A, b)
215 %!shared A, b, n, M1, M2
219 %! A = spdiags ([-2*ones(n,1) 4*ones(n,1) -ones(n,1)], -1:1, n, n);
223 %! M1 = spdiags ([ones(n,1)/(-2) ones(n,1)],-1:0, n, n);
224 %! M2 = spdiags ([4*ones(n,1) -ones(n,1)], 0:1, n, n);
225 %! [x, flag, relres, iter, resvec] = bicgstab (A, b, tol, maxit, M1, M2);
226 %! assert (x, ones (size (b)), 1e-7);
229 %!function y = afun (x, a)
236 %! [x, flag, relres, iter, resvec] = bicgstab (@(x) afun (x, A), b,
237 %! tol, maxit, M1, M2);
238 %! assert (x, ones (size (b)), 1e-7);
243 %! a = sprand (n, n, .1);
244 %! A = a'*a + 100 * eye (n);
246 %! [x, flag, relres, iter, resvec] = bicgstab (A, b, tol, [], diag (diag (A)));
247 %! assert (x, ones (size (b)), 1e-7);