]> Creatis software - CreaPhase.git/blob - octave_packages/m/sparse/bicgstab.m
update packages
[CreaPhase.git] / octave_packages / m / sparse / bicgstab.m
1 ## Copyright (C) 2008-2012 Radek Salac
2 ## Copyright (C) 2012 Carlo de Falco
3 ##
4 ## This file is part of Octave.
5 ##
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.
10 ##
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.
15 ##
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/>.
19
20 ## -*- texinfo -*-
21 ##
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
26 ## method.
27 ##
28 ## @itemize @minus
29 ## @item @var{rtol} is the relative tolerance, if not given or set to
30 ## [] the default value 1e-6 is used.
31 ##
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
34 ## used.
35 ##
36 ## @item @var{x0} the initial guess, if not given or set to [] the
37 ## default value @code{zeros (size (b))} is used.
38 ## @end itemize
39 ##
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}.
42 ##
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}.
47 ##
48 ## If called with more than one output parameter
49 ##
50 ## @itemize @minus
51 ## @item @var{flag} indicates the exit status:
52 ## @itemize @minus
53 ## @item 0: iteration converged to the within the chosen tolerance
54 ##
55 ## @item 1: the maximum number of iterations was reached before convergence
56 ##
57 ## @item 3: the algorithm reached stagnation
58 ## @end itemize
59 ## (the value 2 is unused but skipped for compatibility).
60 ##
61 ## @item @var{relres} is the final value of the relative residual.
62 ##
63 ## @item @var{iter} is the number of iterations performed.
64 ##
65 ## @item @var{resvec} is a vector containing the relative residual at each iteration.
66 ## @end itemize
67 ##
68 ## @seealso{bicg, cgs, gmres, pcg}
69 ##
70 ## @end deftypefn
71
72 function [x, flag, relres, iter, resvec] = bicgstab (A, b, tol, maxit,
73                                                      M1, M2, x0)
74
75   if (nargin >= 2 && nargin <= 7 && isvector (full (b)))
76
77     if (ischar (A))
78       A = str2func (A);
79     elseif (ismatrix (A))
80       Ax  = @(x) A  * x;
81     elseif (isa (A, "function_handle"))
82       Ax  = @(x) feval (A, x);
83     else
84       error (["bicgstab: first argument is expected " ...
85               "to be a function or a square matrix"]);
86     endif
87
88     if (nargin < 3 || isempty (tol))
89       tol = 1e-6;
90     endif
91
92     if (nargin < 4 || isempty (maxit))
93       maxit = min (rows (b), 20);
94     endif
95
96     if (nargin < 5 || isempty (M1))
97       M1m1x = @(x) x;
98     elseif (ischar (M1))
99       M1m1x = str2func (M1);
100     elseif (ismatrix (M1))
101       M1m1x = @(x) M1  \ x;
102     elseif (isa (M1, "function_handle"))
103       M1m1x = @(x) feval (M1, x);
104     else
105       error (["bicgstab: preconditioner is " ...
106               "expected to be a function or matrix"]);
107     endif
108
109     if (nargin < 6 || isempty (M2))
110       M2m1x = @(x) x;
111     elseif (ischar (M2))
112       M2m1x = str2func (M2);
113     elseif (ismatrix (M2))
114       M2m1x = @(x) M2  \ x;
115     elseif (isa (M2, "function_handle"))
116       M2m1x = @(x) feval (M2, x);
117     else
118       error (["bicgstab: preconditioner is "...
119               "expected to be a function or matrix"]);
120     endif
121
122     precon = @(x) M2m1x (M1m1x (x));
123
124     if (nargin < 7 || isempty (x0))
125       x0 = zeros (size (b));
126     endif
127
128     ## specifies initial estimate x0
129     if (nargin < 7)
130       x = zeros (rows (b), 1);
131     else
132       x = x0;
133     endif
134
135     norm_b = norm (b);
136
137     res = b - Ax (x);
138     rr = res;
139
140     ## Vector of the residual norms for each iteration.
141     resvec = norm(res) / norm_b;
142
143     ## Default behaviour we don't reach tolerance tol within maxit iterations.
144     flag = 1;
145
146     for iter = 1:maxit
147       rho_1 = res' * rr;
148
149       if (iter == 1)
150         p = res;
151       else
152         beta = (rho_1 / rho_2) * (alpha / omega);
153         p = res + beta * (p - omega * v);
154       endif
155
156       phat = precon (p);
157
158       v = Ax (phat);
159       alpha = rho_1 / (rr' * v);
160       s = res - alpha * v;
161
162       shat = precon (s);
163
164       t = Ax (shat);
165       omega = (t' * s) / (t' * t);
166       x = x + alpha * phat + omega * shat;
167       res = s - omega * t;
168       rho_2 = rho_1;
169
170       relres = norm (res) / norm_b;
171       resvec = [resvec; relres];
172
173       if (relres <= tol)
174         ## We reach tolerance tol within maxit iterations.
175         flag = 0;
176         break;
177       elseif (resvec(end) == resvec(end - 1))
178         ## The method stagnates.
179         flag = 3;
180         break;
181       endif
182     endfor
183
184     if (nargout < 2)
185       if (flag == 0)
186         printf ("bicgstab converged at iteration %i ", iter);
187         printf ("to a solution with relative residual %e\n", relres);
188       elseif (flag == 3)
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);
194       else
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);
200       endif
201     endif
202
203   else
204     print_usage ();
205   endif
206
207 endfunction
208
209 %!demo
210 %! % Solve system of A*x=b
211 %! A = [5 -1 3;-1 2 -2;3 -2 3]
212 %! b = [7;-1;4]
213 %! [x, flag, relres, iter, resvec] = bicgstab(A, b)
214
215 %!shared A, b, n, M1, M2
216 %!
217 %!test
218 %! n = 100;
219 %! A = spdiags ([-2*ones(n,1) 4*ones(n,1) -ones(n,1)], -1:1, n, n);
220 %! b = sum (A, 2);
221 %! tol = 1e-8;
222 %! maxit = 15;
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);
227 %!
228 %!test
229 %!function y = afun (x, a)
230 %!  y = a * x;
231 %!endfunction
232 %!
233 %! tol = 1e-8;
234 %! maxit = 15;
235 %!
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);
239
240 %!test
241 %! n = 100;
242 %! tol = 1e-8;
243 %! a = sprand (n, n, .1);
244 %! A = a'*a + 100 * eye (n);
245 %! b = sum (A, 2);
246 %! [x, flag, relres, iter, resvec] = bicgstab (A, b, tol, [], diag (diag (A)));
247 %! assert (x, ones (size (b)), 1e-7);