Mercurial > hg > octave-lyh
annotate scripts/sparse/bicg.m @ 14853:72b8b39e12be
doc: Periodic grammarcheck of documentation.
* contrib.txi, diagperm.txi, emacs.txi, install.txi, package.txi, plot.txi,
poly.txi, vectorize.txi, strread.m, textscan.m, graphics_toolkit.m, bicg.m,
bicgstab.m, cgs.m, rand.cc, data.cc: Periodic grammarcheck of documentation.
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Mon, 09 Jul 2012 10:34:43 -0700 |
parents | 3edba8b5f430 |
children | e2de3c8882be |
rev | line source |
---|---|
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
1 ## Copyright (C) 2006 Sylvain Pelissier <sylvain.pelissier@gmail.com> |
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
13931
diff
changeset
|
2 ## Copyright (C) 2012 Carlo de Falco |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
3 ## |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
4 ## This program is free software; you can redistribute it and/or modify |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
5 ## it under the terms of the GNU General Public License as published by |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
6 ## the Free Software Foundation; either version 2 of the License, or |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
7 ## (at your option) any later version. |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
8 ## |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
9 ## This program is distributed in the hope that it will be useful, |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
10 ## but WITHOUT ANY WARRANTY; without even the implied warranty of |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
11 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
12 ## GNU General Public License for more details. |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
13 ## |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
14 ## You should have received a copy of the GNU General Public License |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
15 ## along with this program; If not, see <http://www.gnu.org/licenses/>. |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
16 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
17 ## -*- texinfo -*- |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
18 ## @deftypefn {Function File} {@var{x} =} bicg (@var{A}, @var{b}, @var{rtol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0}) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
19 ## @deftypefnx {Function File} {@var{x} =} bicg (@var{A}, @var{b}, @var{rtol}, @var{maxit}, @var{P}) |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
20 ## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} bicg (@var{A}, @var{b}, @dots{}) |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
21 ## Solve @code{A x = b} using the Bi-conjugate gradient iterative method. |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
22 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
23 ## @itemize @minus |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
24 ## @item @var{rtol} is the relative tolerance, if not given |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
25 ## or set to [] the default value 1e-6 is used. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
26 ## |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
27 ## @item @var{maxit} the maximum number of outer iterations, |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
28 ## if not given or set to [] the default value |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
29 ## @code{min (20, numel (b))} is used. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
30 ## |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
31 ## @item @var{x0} the initial guess, if not given or set to [] |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
32 ## the default value @code{zeros (size (b))} is used. |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
33 ## @end itemize |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
34 ## |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
35 ## @var{A} can be passed as a matrix or as a function handle or |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
36 ## inline function @code{f} such that @code{f(x, "notransp") = A*x} |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
37 ## and @code{f(x, "transp") = A'*x}. |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
38 ## |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
39 ## The preconditioner @var{P} is given as @code{P = M1 * M2}. |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
40 ## Both @var{M1} and @var{M2} can be passed as a matrix or as |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
41 ## a function handle or inline function @code{g} such that |
14359
7277fe922e99
doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
42 ## @code{g(x, "notransp") = M1 \ x} or @code{g(x, "notransp") = M2 \ x} and |
7277fe922e99
doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
43 ## @code{g(x, "transp") = M1' \ x} or @code{g(x, "transp") = M2' \ x}. |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
44 ## |
13931
9de488c6c59c
doc: Spellcheck documentation before 3.6.0 release
Rik <octave@nomad.inbox5.com>
parents:
13929
diff
changeset
|
45 ## If called with more than one output parameter |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
46 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
47 ## @itemize @minus |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
48 ## @item @var{flag} indicates the exit status: |
14853
72b8b39e12be
doc: Periodic grammarcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14609
diff
changeset
|
49 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
50 ## @itemize @minus |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
51 ## @item 0: iteration converged to the within the chosen tolerance |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
52 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
53 ## @item 1: the maximum number of iterations was reached before convergence |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
54 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
55 ## @item 3: the algorithm reached stagnation |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
56 ## @end itemize |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
57 ## (the value 2 is unused but skipped for compatibility). |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
58 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
59 ## @item @var{relres} is the final value of the relative residual. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
60 ## |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
61 ## @item @var{iter} is the number of iterations performed. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
62 ## |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
63 ## @item @var{resvec} is a vector containing the relative residual at each iteration. |
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
64 ## @end itemize |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
65 ## |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
66 ## @seealso{bicgstab, cgs, gmres, pcg} |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
67 ## |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
68 ## @end deftypefn |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
69 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
70 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
71 function [x, flag, res1, k, resvec] = bicg (A, b, tol, maxit, M1, M2, x0) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
72 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
73 if (nargin >= 2 && isvector (full (b))) |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
74 |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
75 if (ischar (A)) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
76 fun = str2func (A); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
77 Ax = @(x) feval (fun, x, "notransp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
78 Atx = @(x) feval (fun, x, "transp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
79 elseif (ismatrix (A)) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
80 Ax = @(x) A * x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
81 Atx = @(x) A' * x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
82 elseif (isa (A, "function_handle")) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
83 Ax = @(x) feval (A, x, "notransp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
84 Atx = @(x) feval (A, x, "transp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
85 else |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
86 error (["bicg: first argument is expected to " ... |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
87 "be a function or a square matrix"]); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
88 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
89 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
90 if (nargin < 3 || isempty (tol)) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
91 tol = 1e-6; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
92 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
93 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
94 if (nargin < 4 || isempty (maxit)) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
95 maxit = min (rows (b), 20); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
96 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
97 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
98 if (nargin < 5 || isempty (M1)) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
99 M1m1x = @(x, ignore) x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
100 M1tm1x = M1m1x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
101 elseif (ischar (M1)) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
102 fun = str2func (M1); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
103 M1m1x = @(x) feval (fun, x, "notransp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
104 M1tm1x = @(x) feval (fun, x, "transp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
105 elseif (ismatrix (M1)) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
106 M1m1x = @(x) M1 \ x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
107 M1tm1x = @(x) M1' \ x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
108 elseif (isa (M1, "function_handle")) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
109 M1m1x = @(x) feval (M1, x, "notransp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
110 M1tm1x = @(x) feval (M1, x, "transp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
111 else |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
112 error (["bicg: preconditioner is expected to " ... |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
113 "be a function or matrix"]); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
114 endif |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
115 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
116 if (nargin < 6 || isempty (M2)) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
117 M2m1x = @(x, ignore) x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
118 M2tm1x = M2m1x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
119 elseif (ischar (M2)) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
120 fun = str2func (M2); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
121 M2m1x = @(x) feval (fun, x, "notransp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
122 M2tm1x = @(x) feval (fun, x, "transp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
123 elseif (ismatrix (M2)) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
124 M2m1x = @(x) M2 \ x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
125 M2tm1x = @(x) M2' \ x; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
126 elseif (isa (M2, "function_handle")) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
127 M2m1x = @(x) feval (M2, x, "notransp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
128 M2tm1x = @(x) feval (M2, x, "transp"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
129 else |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
130 error (["bicg: preconditioner is expected to " ... |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
131 "be a function or matrix"]); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
132 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
133 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
134 Pm1x = @(x) M2m1x (M1m1x (x)); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
135 Ptm1x = @(x) M1tm1x (M2tm1x (x)); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
136 |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
137 if (nargin < 7 || isempty (x0)) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
138 x0 = zeros (size (b)); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
139 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
140 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
141 y = x = x0; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
142 c = b; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
143 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
144 r0 = b - Ax (x); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
145 s0 = c - Atx (y); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
146 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
147 d = Pm1x (r0); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
148 f = Ptm1x (s0); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
149 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
150 bnorm = norm (b); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
151 res0 = Inf; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
152 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
153 if (any (r0 != 0)) |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
154 |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
155 for k = 1:maxit |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
156 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
157 a = (s0' * Pm1x (r0)) ./ (f' * Ax (d)); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
158 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
159 x += a * d; |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
160 y += conj (a) * f; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
161 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
162 r1 = r0 - a * Ax (d); |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
163 s1 = s0 - conj (a) * Atx (f); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
164 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
165 beta = (s1' * Pm1x (r1)) ./ (s0' * Pm1x (r0)); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
166 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
167 d = Pm1x (r1) + beta * d; |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
168 f = Ptm1x (s1) + conj (beta) * f; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
169 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
170 r0 = r1; |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
171 s0 = s1; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
172 |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
173 res1 = norm (b - Ax (x)) / bnorm; |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
174 if (res1 < tol) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
175 flag = 0; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
176 if (nargout < 2) |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
177 printf ("bicg converged at iteration %i ", k); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
178 printf ("to a solution with relative residual %e\n", res1); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
179 endif |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
180 break; |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
181 endif |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
182 |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
183 if (res0 <= res1) |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
184 flag = 3; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
185 printf ("bicg stopped at iteration %i ", k); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
186 printf ("without converging to the desired tolerance %e\n", tol); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
187 printf ("because the method stagnated.\n"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
188 printf ("The iterate returned (number %i) ", k-1); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
189 printf ("has relative residual %e\n", res0); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
190 break |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
191 endif |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
192 res0 = res1; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
193 if (nargout > 4) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
194 resvec(k) = res0; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
195 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
196 endfor |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
197 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
198 if (k == maxit) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
199 flag = 1; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
200 printf ("bicg stopped at iteration %i ", maxit); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
201 printf ("without converging to the desired tolerance %e\n", tol); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
202 printf ("because the maximum number of iterations was reached. "); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
203 printf ("The iterate returned (number %i) has ", maxit); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
204 printf ("relative residual %e\n", res1); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
205 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
206 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
207 else |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
208 flag = 0; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
209 if (nargout < 2) |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
210 printf ("bicg converged after 0 interations\n"); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
211 endif |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
212 endif |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
213 |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
214 else |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
215 print_usage (); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
216 endif |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
217 |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
218 endfunction; |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
219 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
220 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
221 %!test |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
222 %! n = 100; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
223 %! A = spdiags ([-2*ones(n,1) 4*ones(n,1) -ones(n,1)], -1:1, n, n); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
224 %! b = sum (A, 2); |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
225 %! tol = 1e-8; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
226 %! maxit = 15; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
227 %! M1 = spdiags ([ones(n,1)/(-2) ones(n,1)],-1:0, n, n); |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
228 %! M2 = spdiags ([4*ones(n,1) -ones(n,1)], 0:1, n, n); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
229 %! [x, flag, relres, iter, resvec] = bicg (A, b, tol, maxit, M1, M2); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
230 %! assert (x, ones (size (b)), 1e-7); |
13796
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
231 |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
232 %!function y = afun (x, t, a) |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
233 %! switch t |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
234 %! case "notransp" |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
235 %! y = a * x; |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
236 %! case "transp" |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
237 %! y = a' * x; |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
238 %! endswitch |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
239 %!endfunction |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
240 %! |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
241 %!test |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
242 %! n = 100; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
243 %! A = spdiags ([-2*ones(n,1) 4*ones(n,1) -ones(n,1)], -1:1, n, n); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
244 %! b = sum (A, 2); |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
245 %! tol = 1e-8; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
246 %! maxit = 15; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
247 %! M1 = spdiags ([ones(n,1)/(-2) ones(n,1)],-1:0, n, n); |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
248 %! M2 = spdiags ([4*ones(n,1) -ones(n,1)], 0:1, n, n); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
249 %! |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
250 %! [x, flag, relres, iter, resvec] = bicg (@(x, t) afun (x, t, A), |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
251 %! b, tol, maxit, M1, M2); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
252 %! assert (x, ones (size (b)), 1e-7); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
253 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
254 %!test |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
255 %! n = 100; |
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
256 %! tol = 1e-8; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
257 %! a = sprand (n, n, .1); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
258 %! A = a' * a + 100 * eye (n); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
259 %! b = sum (A, 2); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
260 %! [x, flag, relres, iter, resvec] = bicg (A, b, tol, [], diag (diag (A))); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
261 %! assert (x, ones (size (b)), 1e-7); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14359
diff
changeset
|
262 |