Mercurial > octave
annotate scripts/sparse/bicg.m @ 30893:e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
* accumarray.m, accumdim.m, quadl.m, quadv.m, randi.m, structfun.m,
__is_function__.m, uigetfile.m, uimenu.m, uiputfile.m, doc_cache_create.m,
colorspace_conversion_input_check.m, imageIO.m, argnames.m, vectorize.m,
vectorize.m, normest1.m, inputname.m, nthargout.m, display_info_file.m,
decic.m, ode15i.m, ode15s.m, ode23.m, ode23s.m, ode45.m, odeset.m,
check_default_input.m, integrate_adaptive.m, ode_event_handler.m,
runge_kutta_23.m, runge_kutta_23s.m, runge_kutta_45_dorpri.m,
runge_kutta_interpolate.m, starting_stepsize.m, __all_opts__.m, fminbnd.m,
fminsearch.m, fminunc.m, fsolve.m, fzero.m, sqp.m, fplot.m, plotyy.m,
__bar__.m, __ezplot__.m, flat_entry.html, profexport.m, movfun.m, bicg.m,
bicgstab.m, cgs.m, eigs.m, gmres.m, pcg.m, __alltohandles__.m, __sprand__.m,
qmr.m, tfqmr.m, dump_demos.m:
Replace "func", "fun", "fn" in documentation and variable names with "fcn".
author | Rik <rik@octave.org> |
---|---|
date | Mon, 04 Apr 2022 18:14:56 -0700 |
parents | 796f54d4ddbf |
children | 597f3ee61a48 |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 ######################################################################## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 ## |
30564
796f54d4ddbf
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
29359
diff
changeset
|
3 ## Copyright (C) 2006-2022 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
27356
diff
changeset
|
4 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 ## See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 ## distribution or <https://octave.org/copyright/>. |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
7 ## |
17795
0a8c35ae5ce1
maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
8 ## This file is part of Octave. |
0a8c35ae5ce1
maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
9 ## |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23084
diff
changeset
|
10 ## Octave is free software: you can redistribute it and/or modify it |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
11 ## under the terms of the GNU General Public License as published by |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23084
diff
changeset
|
12 ## the Free Software Foundation, either version 3 of the License, or |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
13 ## (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
|
14 ## |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
15 ## Octave is distributed in the hope that it will be useful, but |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
17 ## 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
|
18 ## 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
|
19 ## |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
20 ## You should have received a copy of the GNU General Public License |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
21 ## along with Octave; see the file COPYING. If not, see |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23084
diff
changeset
|
22 ## <https://www.gnu.org/licenses/>. |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
23 ## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 ######################################################################## |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
25 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
26 ## -*- texinfo -*- |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
27 ## @deftypefn {} {@var{x} =} bicg (@var{A}, @var{b}) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
28 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
29 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}, @var{maxit}) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
30 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M}) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
31 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M1}, @var{M2}) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
32 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M}, [], @var{x0}) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
33 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0}) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
34 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M}, [], @var{x0}, @dots{}) |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
35 ## @deftypefnx {} {@var{x} =} bicg (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0}, @dots{}) |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20711
diff
changeset
|
36 ## @deftypefnx {} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} bicg (@var{A}, @var{b}, @dots{}) |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
37 ## Solve the linear system of equations @w{@code{@var{A} * @var{x} = @var{b}}} |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
38 ## by means of the Bi-Conjugate Gradient iterative method. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
39 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
40 ## The input arguments are: |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
41 ## |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
42 ## @itemize |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
43 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
44 ## @item @var{A} is the matrix of the linear system and it must be square. |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
45 ## @var{A} can be passed as a matrix, function handle, or inline function |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
46 ## @code{Afcn} such that @w{@code{Afcn (x, "notransp") = A * x}} and |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
47 ## @w{@code{Afcn (x, "transp") = A' * x}}. Additional parameters to |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
48 ## @code{Afcn} may be passed after @var{x0}. |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
49 ## |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
50 ## @item @var{b} is the right-hand side vector. It must be a column vector |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
51 ## with the same number of rows as @var{A}. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
52 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
53 ## @item |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
54 ## @var{tol} is the required relative tolerance for the residual error, |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
55 ## @w{@code{@var{b} - @var{A} * @var{x}}}. The iteration stops if |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
56 ## @w{@code{norm (@var{b} - @var{A} * @var{x})} @leq{} |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
57 ## @w{@code{@var{tol} * norm (@var{b})}}}. |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
58 ## If @var{tol} is omitted or empty, then a tolerance of 1e-6 is used. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
59 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
60 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
61 ## @var{maxit} is the maximum allowed number of iterations; if @var{maxit} |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
62 ## is omitted or empty then a value of 20 is used. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
63 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
64 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
65 ## @var{M1}, @var{M2} are the preconditioners. The preconditioner @var{M} is |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
66 ## given as @code{@var{M} = @var{M1} * @var{M2}}. Both @var{M1} and @var{M2} |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
67 ## can be passed as a matrix or as a function handle or inline function |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
68 ## @code{g} such that @w{@code{g (@var{x}, "notransp") = @var{M1} \ @var{x}}} |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
69 ## or @w{@code{g (@var{x}, "notransp") = @var{M2} \ @var{x}}} and |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
70 ## @w{@code{g (@var{x}, "transp") = @var{M1}' \ @var{x}}} or |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
71 ## @w{@code{g (@var{x}, "transp") = @var{M2}' \ @var{x}}}. |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
72 ## If @var{M1} is omitted or empty, then preconditioning is not applied. |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
73 ## The preconditioned system is theoretically equivalent to applying the |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
74 ## @code{bicg} method to the linear system |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
75 ## @code{inv (@var{M1}) * A * inv (@var{M2}) * @var{y} = inv |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
76 ## (@var{M1}) * @var{b}} and |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
77 ## @code{inv (@var{M2'}) * A' * inv (@var{M1'}) * @var{z} = |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
78 ## inv (@var{M2'}) * @var{b}} and then setting |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
79 ## @code{@var{x} = inv (@var{M2}) * @var{y}}. |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
80 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
81 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
82 ## @var{x0} is the initial guess. If @var{x0} is omitted or empty then the |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
83 ## function sets @var{x0} to a zero vector by default. |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
84 ## @end itemize |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
85 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
86 ## Any arguments which follow @var{x0} are treated as parameters, and passed in |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
87 ## an appropriate manner to any of the functions (@var{Afcn} or @var{Mfcn}) or |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
88 ## that have been given to @code{bicg}. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
89 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
90 ## The output parameters are: |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
91 ## |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
92 ## @itemize |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
93 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
94 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
95 ## @var{x} is the computed approximation to the solution of |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
96 ## @w{@code{@var{A} * @var{x} = @var{b}}}. If the algorithm did not converge, |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
97 ## then @var{x} is the iteration which has the minimum residual. |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
98 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
99 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
100 ## @var{flag} indicates the exit status: |
14853
72b8b39e12be
doc: Periodic grammarcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14609
diff
changeset
|
101 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
102 ## @itemize |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
103 ## @item 0: The algorithm converged to within the prescribed tolerance. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
104 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
105 ## @item 1: The algorithm did not converge and it reached the maximum number of |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
106 ## iterations. |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
107 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
108 ## @item 2: The preconditioner matrix is singular. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
109 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
110 ## @item 3: The algorithm stagnated, i.e., the absolute value of the |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
111 ## difference between the current iteration @var{x} and the previous is less |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
112 ## than @code{eps * norm (@var{x},2)}. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
113 ## |
25200
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
114 ## @item 4: The algorithm could not continue because intermediate values |
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
115 ## became too small or too large for reliable computation. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
116 ## @end itemize |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
117 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
118 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
119 ## @var{relres} is the ratio of the final residual to its initial value, |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
120 ## measured in the Euclidean norm. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
121 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
122 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
123 ## @var{iter} is the iteration which @var{x} is computed. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
124 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
125 ## @item |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
126 ## @var{resvec} is a vector containing the residual at each iteration. |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
127 ## The total number of iterations performed is given by |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
128 ## @code{length (@var{resvec}) - 1}. |
13091
e5aaba072d2b
maint: style fixes in iterative linear solvers
Carlo de Falco <kingcrimson@tiscali.it>
parents:
13073
diff
changeset
|
129 ## @end itemize |
16816
12005245b645
doc: Periodic grammarcheck of documentation.
Rik <rik@octave.org>
parents:
16768
diff
changeset
|
130 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
131 ## Consider a trivial problem with a tridiagonal matrix |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
132 ## |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
133 ## @example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
134 ## @group |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
135 ## n = 20; |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
136 ## A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n)) + ... |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
137 ## toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ... |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
138 ## sparse (1, 2, 1, 1, n) * n / 2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
139 ## b = A * ones (n, 1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
140 ## restart = 5; |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
141 ## [M1, M2] = ilu (A); # in this tridiag case, it corresponds to lu (A) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
142 ## M = M1 * M2; |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
143 ## Afcn = @@(x, string) strcmp (string, "notransp") * (A * x) + ... |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
144 ## strcmp (string, "transp") * (A' * x); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
145 ## Mfcn = @@(x, string) strcmp (string, "notransp") * (M \ x) + ... |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
146 ## strcmp (string, "transp") * (M' \ x); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
147 ## M1fcn = @@(x, string) strcmp (string, "notransp") * (M1 \ x) + ... |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
148 ## strcmp (string, "transp") * (M1' \ x); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
149 ## M2fcn = @@(x, string) strcmp (string, "notransp") * (M2 \ x) + ... |
25579
07c2c42f457e
doc: Miscellaneous documentation fixes all over the manual (bug #54288).
Rik <rik@octave.org>
parents:
25200
diff
changeset
|
150 ## strcmp (string, "transp") * (M2' \ x); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
151 ## @end group |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
152 ## @end example |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
153 ## |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
154 ## @sc{Example 1:} simplest usage of @code{bicg} |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
155 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
156 ## @example |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
157 ## x = bicg (A, b) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
158 ## @end example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
159 ## |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
160 ## @sc{Example 2:} @code{bicg} with a function that computes |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
161 ## @code{@var{A}*@var{x}} and @code{@var{A'}*@var{x}} |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
162 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
163 ## @example |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
164 ## x = bicg (Afcn, b, [], n) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
165 ## @end example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
166 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
167 ## @sc{Example 3:} @code{bicg} with a preconditioner matrix @var{M} |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
168 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
169 ## @example |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
170 ## x = bicg (A, b, 1e-6, n, M) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
171 ## @end example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
172 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
173 ## @sc{Example 4:} @code{bicg} with a function as preconditioner |
13929
9cae456085c2
Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents:
13796
diff
changeset
|
174 ## |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
175 ## @example |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
176 ## x = bicg (Afcn, b, 1e-6, n, Mfcn) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
177 ## @end example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
178 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
179 ## @sc{Example 5:} @code{bicg} with preconditioner matrices @var{M1} |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
180 ## and @var{M2} |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
181 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
182 ## @example |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
183 ## x = bicg (A, b, 1e-6, n, M1, M2) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
184 ## @end example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
185 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
186 ## @sc{Example 6:} @code{bicg} with functions as preconditioners |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
187 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
188 ## @example |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
189 ## x = bicg (Afcn, b, 1e-6, n, M1fcn, M2fcn) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
190 ## @end example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
191 ## |
24992
3d5a39077708
avoid some warnings from old versions of makeinfo (bug #53479)
Guillaume Flandin
parents:
24985
diff
changeset
|
192 ## @sc{Example 7:} @code{bicg} with as input a function requiring an argument |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
193 ## |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
194 ## @example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
195 ## @group |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
196 ## function y = Ap (A, x, string, z) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
197 ## ## compute A^z * x or (A^z)' * x |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
198 ## y = x; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
199 ## if (strcmp (string, "notransp")) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
200 ## for i = 1:z |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
201 ## y = A * y; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
202 ## endfor |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
203 ## elseif (strcmp (string, "transp")) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
204 ## for i = 1:z |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
205 ## y = A' * y; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
206 ## endfor |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
207 ## endif |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
208 ## endfunction |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
209 ## |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
210 ## Apfcn = @@(x, string, p) Ap (A, x, string, p); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
211 ## x = bicg (Apfcn, b, [], [], [], [], [], 2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
212 ## @end group |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
213 ## @end example |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
214 ## |
27984
b09432b20a84
maint: Remove special cases of old version control keywords in code base.
Rik <rik@octave.org>
parents:
27923
diff
changeset
|
215 ## Reference: |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
216 ## |
27984
b09432b20a84
maint: Remove special cases of old version control keywords in code base.
Rik <rik@octave.org>
parents:
27923
diff
changeset
|
217 ## @nospell{Y. Saad}, @cite{Iterative Methods for Sparse Linear Systems}, |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
218 ## Second edition, 2003, SIAM. |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
219 ## |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
220 ## @seealso{bicgstab, cgs, gmres, pcg, qmr, tfqmr} |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
221 ## @end deftypefn |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
222 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
223 function [x_min, flag, relres, iter_min, resvec] = ... |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
224 bicg (A, b, tol = [], maxit = [], M1 = [], M2 = [], x0 = [], varargin) |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
225 |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
226 [Afcn, M1fcn, M2fcn] = __alltohandles__ (A, b, M1, M2, "bicg"); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
227 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
228 [tol, maxit, x0] = __default__input__ ({1e-06, min(rows(b), 20), ... |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
229 zeros(rows (b),1)}, tol, maxit, x0); |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
230 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
231 if (columns (b) == 2) |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
232 c = b(:,2); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
233 b = b(:,1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
234 else |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
235 c = b; |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
236 endif |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
237 norm_b = norm (b, 2); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
238 |
28912
0de38a6ef693
maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents:
27984
diff
changeset
|
239 if (norm_b == 0) # the only (only iff det (A) == 0) solution is x = 0 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
240 if (nargout < 2) |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
241 printf ("The right hand side vector is all zero so bicg\n") |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
242 printf ("returned an all zero solution without iterating.\n") |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
243 endif |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
244 x_min = zeros (numel (b), 1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
245 flag = 0; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
246 relres = 0; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
247 iter_min = 0; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
248 resvec = 0; |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
249 return; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
250 endif |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
251 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
252 x = x_min = x_pr = x0; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
253 iter = iter_min = 0; |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
254 flag = 1; # Default flag is "maximum number of iterations reached" |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
255 resvec = zeros (maxit + 1, 1); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
256 r0 = b - Afcn (x, "notransp", varargin{:}); # Residual of the system |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
257 s0 = c - Afcn (x, "transp", varargin{:}); # Res. of the "dual system" |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
258 resvec(1) = norm (r0, 2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
259 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
260 try |
28928
ae7ce8358953
maint: Add semicolon to end of all warning() and error() invocations.
Rik <rik@octave.org>
parents:
28917
diff
changeset
|
261 warning ("error", "Octave:singular-matrix", "local"); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
262 prec_r0 = M1fcn (r0, "notransp", varargin{:}); # r0 preconditioned |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
263 prec_s0 = s0; |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
264 prec_r0 = M2fcn (prec_r0, "notransp", varargin{:}); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
265 prec_s0 = M2fcn (prec_s0, "transp", varargin{:}); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
266 prec_s0 = M1fcn (prec_s0, "transp", varargin{:}); # s0 preconditioned |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
267 p = prec_r0; # Direction of the system |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
268 q = prec_s0; # Direction of the "dual system" |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
269 catch |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
270 flag = 2; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
271 end_try_catch |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
272 |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
273 while ((flag != 2) && (iter < maxit) && (resvec(iter+1) >= norm_b * tol)) |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
274 v = Afcn (p, "notransp", varargin{:}); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
275 prod_qv = q' * v; |
25200
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
276 alpha = (s0' * prec_r0); |
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
277 if (abs (prod_qv) <= eps * abs (alpha)) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
278 flag = 4; |
28917
72d57dbcc305
maint: Add semicolon after break and return keywords.
Rik <rik@octave.org>
parents:
28912
diff
changeset
|
279 break; |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
280 endif |
25200
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
281 alpha ./= prod_qv; |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
282 x += alpha * p; |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
283 prod_rs = (s0' * prec_r0); # Product between r0 and s0 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
284 r0 -= alpha * v; |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
285 s0 -= conj (alpha) * Afcn (q, "transp", varargin{:}); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
286 prec_r0 = M1fcn (r0, "notransp", varargin{:}); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
287 prec_s0 = s0; |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
288 prec_r0 = M2fcn (prec_r0, "notransp", varargin{:}); |
25200
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
289 beta = s0' * prec_r0; |
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
290 if (abs (prod_rs) <= abs (beta)) |
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
291 flag = 4; |
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
292 break; |
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
293 endif |
b3ee0179d7b0
bicg.m: Check for denominators much smaller than numerators bicg.m (bug #53589).
Marco Caliari <marco.caliari@univr.it>
parents:
25166
diff
changeset
|
294 beta ./= prod_rs; |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
295 prec_s0 = M2fcn (prec_s0, "transp", varargin{:}); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
296 prec_s0 = M1fcn (prec_s0, "transp", varargin{:}); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
297 iter += 1; |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
298 resvec(iter+1) = norm (r0); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
299 if (resvec(iter+1) <= resvec(iter_min+1)) |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
300 x_min = x; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
301 iter_min = iter; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
302 endif |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
303 if (norm (x - x_pr) <= norm (x) * eps) |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
304 flag = 3; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
305 break; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
306 endif |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
307 p = prec_r0 + beta*p; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
308 q = prec_s0 + conj (beta) * q; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
309 endwhile |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
310 resvec = resvec(1:iter+1,1); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
311 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
312 if (flag == 2) |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
313 relres = 1; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
314 else |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
315 relres = resvec(iter_min+1) / norm_b; |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
316 endif |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
317 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
318 if ((flag == 1) && (relres <= tol)) |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
319 flag = 0; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
320 endif |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
321 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
322 if (nargout < 2) |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
323 switch (flag) |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
324 case 0 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
325 printf ("bicg converged at iteration %i ", iter_min); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
326 printf ("to a solution with relative residual %e\n", relres); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
327 case 1 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
328 printf ("bicg stopped at iteration %i ", iter); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
329 printf ("without converging to the desired tolerance %e\n", tol); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
330 printf ("because the maximum number of iterations was reached. "); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
331 printf ("The iterate returned (number %i) has ", iter_min); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
332 printf ("relative residual %e\n", relres); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
333 case 2 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
334 printf ("bicg stopped at iteration %i ", iter); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
335 printf ("without converging to the desired tolerance %e\n", tol); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
336 printf ("because the preconditioner matrix is singular.\n"); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
337 printf ("The iterate returned (number %i) ", iter_min); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
338 printf ("has relative residual %e\n", relres); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
339 case 3 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
340 printf ("bicg stopped at iteration %i ", iter); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
341 printf ("without converging to the desired tolerance %e\n", tol); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
342 printf ("because the method stagnated.\n"); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
343 printf ("The iterate returned (number %i) ", iter_min); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
344 printf ("has relative residual %e\n", relres); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
345 case 4 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
346 printf ("bicg stopped at iteration %i ", iter); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
347 printf ("without converging to the desired tolerance %e\n", tol); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
348 printf ("because the method can't continue.\n"); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
349 printf ("The iterate returned (number %i) ", iter_min); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
350 printf ("has relative residual %e\n", relres); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
351 endswitch |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
352 endif |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
353 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
354 endfunction |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
355 |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
356 |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
357 %!demo |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
358 %! ## simplest use case |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
359 %! n = 20; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
360 %! A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n)) + ... |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
361 %! toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ... |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
362 %! sparse (1, 2, 1, 1, n) * n / 2); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
363 %! b = A * ones (n, 1); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
364 %! [M1, M2] = ilu (A + 0.1 * eye (n)); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
365 %! M = M1 * M2; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
366 %! x = bicg (A, b, [], n); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
367 %! function y = Ap (A, x, string, z) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
368 %! ## compute A^z * x or (A^z)' * x |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
369 %! y = x; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
370 %! if (strcmp (string, "notransp")) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
371 %! for i = 1:z |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
372 %! y = A * y; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
373 %! endfor |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
374 %! elseif (strcmp (string, "transp")) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
375 %! for i = 1:z |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
376 %! y = A' * y; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
377 %! endfor |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
378 %! endif |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
379 %! endfunction |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
380 %! |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
381 %! Afcn = @(x, string) Ap (A, x, string, 1); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
382 %! x = bicg (Afcn, b, [], n); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
383 %! x = bicg (A, b, 1e-6, n, M); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
384 %! x = bicg (A, b, 1e-6, n, M1, M2); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
385 %! function y = Mfcn (M, x, string) |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
386 %! if (strcmp (string, "notransp")) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
387 %! y = M \ x; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
388 %! else |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
389 %! y = M' \ x; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
390 %! endif |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
391 %! endfunction |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
392 %! |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
393 %! M1fcn = @(x, string) Mfcn (M, x, string); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
394 %! x = bicg (Afcn, b, 1e-6, n, M1fcn); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
395 %! M1fcn = @(x, string) Mfcn (M1, x, string); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
396 %! M2fcn = @(x, string) Mfcn (M2, x, string); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
397 %! x = bicg (Afcn, b, 1e-6, n, M1fcn, M2fcn); |
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
398 %! Afcn = @(x, string, p) Ap (A, x, string, p); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
399 %! ## Solution of A^2 * x = b |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
400 %! x = bicg (Afcn, b, [], 2*n, [], [], [], 2); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
401 |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
402 %!test |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
403 %! ## Check that all type of inputs work |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
404 %! A = sparse (toeplitz ([2, 1, 0, 0, 0], [2, -1, 0, 0, 0])); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
405 %! b = A * ones (5, 1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
406 %! M1 = diag (sqrt (diag (A))); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
407 %! M2 = M1; |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
408 %! Afcn = @(z, string) strcmp (string, "notransp") * (A * z) + ... |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
409 %! strcmp (string, "transp") * (A' * z); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
410 %! M1_fcn = @(z, string) strcmp (string,"notransp") * (M1 \ z) + ... |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
411 %! strcmp (string, "transp") * (M1' \ z); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
412 %! M2_fcn = @(z, string) strcmp (string, "notransp") * (M2 \ z) + ... |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
413 %! strcmp (string, "transp") * (M2' \ z); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
414 %! [x, flag] = bicg (A, b); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
415 %! assert (flag, 0); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
416 %! [x, flag] = bicg (A, b, [], [], M1, M2); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
417 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
418 %! [x, flag] = bicg (A, b, [], [], M1_fcn, M2_fcn); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
419 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
420 %! [x, flag] = bicg (A, b,[], [], M1_fcn, M2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
421 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
422 %! [x, flag] = bicg (A, b,[], [], M1, M2_fcn); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
423 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
424 %! [x, flag] = bicg (Afcn, b); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
425 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
426 %! [x, flag] = bicg (Afcn, b,[], [], M1, M2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
427 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
428 %! [x, flag] = bicg (Afcn, b,[], [], M1_fcn, M2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
429 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
430 %! [x, flag] = bicg (Afcn, b,[], [], M1, M2_fcn); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
431 %! assert (flag, 0); |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
432 %! [x, flag] = bicg (Afcn, b,[], [], M1_fcn, M2_fcn); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
433 %! assert (flag, 0); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
434 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
435 %!test |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
436 %! n = 100; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
437 %! 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
|
438 %! b = sum (A, 2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
439 %! tol = 1e-8; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
440 %! maxit = 15; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
441 %! 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
|
442 %! M2 = spdiags ([4*ones(n,1) -ones(n,1)], 0:1, n, n); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
443 %! [x, flag, relres, iter, resvec] = bicg (A, b, tol, maxit, M1, M2); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
444 %! assert (norm (b - A*x) / norm (b), 0, tol); |
13796
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
445 |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
446 %!function y = afcn (x, t, a) |
16933
e39f00a32dc7
maint: Use parentheses around condition for switch(),while(),if() statements.
Rik <rik@octave.org>
parents:
16816
diff
changeset
|
447 %! switch (t) |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
448 %! case "notransp" |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
449 %! y = a * x; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
450 %! case "transp" |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
16933
diff
changeset
|
451 %! y = a' * x; |
13796
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
452 %! endswitch |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
453 %!endfunction |
53674ceb9133
maint: fix function definition in test
John W. Eaton <jwe@octave.org>
parents:
13141
diff
changeset
|
454 %! |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
455 %!test |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
456 %! n = 100; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
457 %! 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
|
458 %! b = sum (A, 2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
459 %! tol = 1e-8; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
460 %! maxit = 15; |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
461 %! 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
|
462 %! 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
|
463 %! |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
464 %! [x, flag, relres, iter, resvec] = bicg (@(x, t) afcn (x, t, A), |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
465 %! b, tol, maxit, M1, M2); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
466 %! 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
|
467 |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
468 %!test |
13141
e81ddf9cacd5
maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
13091
diff
changeset
|
469 %! n = 100; |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
470 %! tol = 1e-8; |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
471 %! a = sprand (n, n, .1); |
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
472 %! 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
|
473 %! b = sum (A, 2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
474 %! [x, flag, relres, iter, resvec] = bicg (A, b, tol, [], diag (diag (A))); |
13021
d55d396a9a55
Add an implementation of the biconjugate gradient iterative method
Carlo de Falco <kingcrimson@tiscali.it>
parents:
diff
changeset
|
475 %! assert (x, ones (size (b)), 1e-7); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
476 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
477 %!test |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
478 %! ## Check that if the preconditioner is singular, the method doesn't work |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
479 %! A = sparse (toeplitz ([2, 1, 0, 0, 0], [2, -1, 0, 0, 0])); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
480 %! b = ones (5,1); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
481 %! M = ones (5); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
482 %! [x, flag] = bicg (A, b, [], [], M); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
483 %! assert (flag, 2); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
484 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
485 %!test |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
486 %! ## If A singular, the algorithm doesn't work due to division by zero |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
487 %! A = ones (5); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
488 %! b = [1:5]'; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
489 %! [x, flag] = bicg (A, b); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
490 %! assert (flag, 4); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
491 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
492 %!test |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
493 %! ## test for a complex linear system |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
494 %! A = sparse (toeplitz ([2, 1, 0, 0, 0], [2, -1, 0, 0, 0]) + ... |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
495 %! 1i * toeplitz ([2, 1, 0, 0, 0], [2, -1, 0, 0, 0])); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
496 %! b = sum (A, 2); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
497 %! [x, flag] = bicg (A, b); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
498 %! assert (flag, 0); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
499 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
500 %!test |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
501 %! A = single (1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
502 %! b = 1; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
503 %! [x, flag] = bicg (A, b); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
504 %! assert (class (x), "single"); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
505 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
506 %!test |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
507 %! A = 1; |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
508 %! b = single (1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
509 %! [x, flag] = bicg (A, b); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
510 %! assert (class (x), "single"); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
511 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
512 %!test |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
513 %! A = single (1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
514 %! b = single (1); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
515 %! [x, flag] = bicg (A, b); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
516 %! assert (class (x), "single"); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
517 |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
518 %!test |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
519 %!function y = Afcn (x, trans) |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
520 %! A = sparse (toeplitz ([2, 1, 0, 0], [2, -1, 0, 0])); |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
521 %! if (strcmp (trans, "notransp")) |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
522 %! y = A * x; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
523 %! else |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
524 %! y = A' * x; |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
525 %! endif |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
526 %!endfunction |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
527 %! |
30893
e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents:
30564
diff
changeset
|
528 %! [x, flag] = bicg ("Afcn", [1; 2; 2; 3]); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
529 %! assert (x, ones (4, 1), 1e-6); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
530 |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
531 %!test |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
532 %! ## unpreconditioned residual |
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
533 %! A = sparse (toeplitz ([2, 1, 0, 0, 0], [2, -1, 0, 0, 0])); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
534 %! b = sum (A, 2); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
535 %! M = magic (5); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
536 %! [x, flag, relres] = bicg (A, b, [], 2, M); |
26160
b9debf4436aa
bicg.m: Relax BIST test by 1eps to pass on some systems (bug #55132).
Rik <rik@octave.org>
parents:
25579
diff
changeset
|
537 %! assert (norm (b - A * x) / norm (b), 0, relres + eps); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
538 |
24863
009ed89dd0fa
test: make bicg and pcg tests conditional on CHOLMOD and UMFPACK
Mike Miller <mtmiller@octave.org>
parents:
24860
diff
changeset
|
539 ## Preconditioned technique |
009ed89dd0fa
test: make bicg and pcg tests conditional on CHOLMOD and UMFPACK
Mike Miller <mtmiller@octave.org>
parents:
24860
diff
changeset
|
540 %!testif HAVE_UMFPACK |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
541 %! A = sparse (toeplitz ([2, 1, 0, 0, 0], [2, -1, 0, 0, 0])); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
542 %! b = sum (A, 2); |
25159
b8279dd83664
Turn off warning about sparse lu factorization in bicg BIST test (bug #53390).
Marco Caliari <marco.caliari@univr.it>
parents:
25152
diff
changeset
|
543 %! warning ("off", "Octave:lu:sparse_input", "local"); |
24860
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
544 %! [M1, M2] = lu (A + eye (5)); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
545 %! [x, flag] = bicg (A, b, [], 1, M1, M2); |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
546 %! ## b has two columns! |
6266e321ef22
Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
24534
diff
changeset
|
547 %! [y, flag] = bicg (M1 \ A / M2, [M1 \ b, M2' \ b], [], 1); |
25152
ad37f04a2eb7
bicg.m: Overhaul GSOC-improved code to conform to Octave conventions.
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
548 %! assert (x, M2 \ y, 8 * eps); |