annotate scripts/sparse/tfqmr.m @ 31225:3eab70385569

sparse-xpow.cc: Use faster multiplication technique, this time for complex
author Arun Giridhar <arungiridhar@gmail.com>
date Sun, 11 Sep 2022 13:53:38 -0400
parents e1788b1a315f
children 597f3ee61a48
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
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) 2016-2022 The Octave Project Developers
27918
b442ec6dda5c use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents: 26376
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/>.
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
7 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
8 ## This file is part of Octave.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
9 ##
24892
2e94172bc1c1 maint: update GPL license text and license header blocks
Mike Miller <mtmiller@octave.org>
parents: 24860
diff changeset
10 ## Octave is free software: you can redistribute it and/or modify it
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
11 ## under the terms of the GNU General Public License as published by
24892
2e94172bc1c1 maint: update GPL license text and license header blocks
Mike Miller <mtmiller@octave.org>
parents: 24860
diff changeset
12 ## the Free Software Foundation, either version 3 of the License, or
2e94172bc1c1 maint: update GPL license text and license header blocks
Mike Miller <mtmiller@octave.org>
parents: 24860
diff changeset
13 ## (at your option) any later version.
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
14 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
15 ## Octave is distributed in the hope that it will be useful, but
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
18 ## GNU General Public License for more details.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
19 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
20 ## You should have received a copy of the GNU General Public License
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
21 ## along with Octave; see the file COPYING. If not, see
24892
2e94172bc1c1 maint: update GPL license text and license header blocks
Mike Miller <mtmiller@octave.org>
parents: 24860
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 ########################################################################
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
25
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
26 ## -*- texinfo -*-
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
27 ## @deftypefn {} {@var{x} =} tfqmr (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0}, @dots{})
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
28 ## @deftypefnx {} {@var{x} =} tfqmr (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M}, [], @var{x0}, @dots{})
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
29 ## @deftypefnx {} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} tfqmr (@var{A}, @var{b}, @dots{})
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
30 ## Solve @code{A x = b} using the Transpose-Tree qmr method, based on the cgs.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
31 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
32 ## The input parameters are:
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
33 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
34 ## @itemize @minus
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
35 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
36 ## @item @var{A} is the matrix of the linear system and it must be square.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
37 ## @var{A} can be passed as a matrix, function handle, or inline
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
38 ## function @code{Afcn} such that @code{Afcn(x) = A * x}. Additional
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
39 ## parameters to @code{Afcn} are passed after @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:
diff changeset
40 ##
24985
d85b2485af9e doc: grammarcheck m-files ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24951
diff changeset
41 ## @item @var{b} is the right hand side vector. It must be a column vector
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
42 ## with the same number of rows as @var{A}.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
43 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
44 ## @item @var{tol} is the relative tolerance, if not given or set to [] the
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
45 ## default value 1e-6 is used.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
46 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
47 ## @item @var{maxit} the maximum number of outer iterations, if not given or
24985
d85b2485af9e doc: grammarcheck m-files ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24951
diff changeset
48 ## set to [] the default value @code{min (20, numel (b))} is used. To be
25003
2365c2661b3c doc: Spellcheck documentation ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24992
diff changeset
49 ## compatible, since the method as different behaviors in the iteration
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
50 ## number is odd or even, is considered as iteration in @code{tfqmr} the
24985
d85b2485af9e doc: grammarcheck m-files ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24951
diff changeset
51 ## entire odd-even cycle. That is, to make an entire iteration, the algorithm
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
52 ## performs two sub-iterations: the odd one and the even one.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
53 ##
24985
d85b2485af9e doc: grammarcheck m-files ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24951
diff changeset
54 ## @item @var{M1}, @var{M2} are the preconditioners. The preconditioner
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
55 ## @var{M} is given as @code{M = M1 * M2}.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
56 ## Both @var{M1} and @var{M2} can be passed as a matrix or as a function
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
57 ## handle or inline function @code{g} such that @code{g(x) = M1 \ x} or
24992
3d5a39077708 avoid some warnings from old versions of makeinfo (bug #53479)
Guillaume Flandin
parents: 24985
diff changeset
58 ## @code{g(x) = M2 \ x}.
25003
2365c2661b3c doc: Spellcheck documentation ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24992
diff changeset
59 ## The technique used is the right-preconditioning, i.e., it is solved
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
60 ## @code{A*inv(M)*y = b} and then @code{x = inv(M)*y}, instead of
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
61 ## @code{A x = b}.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
62 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
63 ## @item @var{x0} the initial guess, if not given or set to [] the default
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
64 ## value @code{zeros (size (b))} is used.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
65 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
66 ## @end itemize
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
67 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
68 ## The arguments which follow @var{x0} are treated as parameters, and passed in
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
69 ## a proper way to any of the functions (@var{A} or @var{M}) which are passed
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
70 ## to @code{tfqmr}.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
71 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
72 ## The output parameters are:
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
73 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
74 ## @itemize @minus
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
75 ##
24985
d85b2485af9e doc: grammarcheck m-files ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24951
diff changeset
76 ## @item @var{x} is the approximation computed. If the method doesn't
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
77 ## converge then it is the iterated with the minimum residual.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
78 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
79 ## @item @var{flag} indicates the exit status:
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
80 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
81 ## @itemize @minus
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
82 ## @item 0: iteration converged to the within the chosen tolerance
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
83 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
84 ## @item 1: the maximum number of iterations was reached before convergence
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
85 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
86 ## @item 2: the preconditioner matrix is singular
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
87 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
88 ## @item 3: the algorithm reached stagnation
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
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:
diff changeset
90 ## @item 4: the algorithm can't continue due to a 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:
diff changeset
91 ## @end itemize
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
92 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
93 ## @item @var{relres} is the relative residual obtained as
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
94 ## @code{(@var{A}*@var{x}-@var{b}) / @code{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:
diff changeset
95 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
96 ## @item @var{iter} is the iteration which @var{x} is
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
97 ## computed.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
98 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
99 ## @item @var{resvec} is a vector containing the residual at each iteration
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
100 ## (including @code{norm (b - A x0)}).
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
101 ## Doing @code{length (@var{resvec}) - 1} is possible to see the
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
102 ## total number of iterations performed.
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
103 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
104 ## @end itemize
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
105 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
106 ## Let us consider a trivial problem with a tridiagonal matrix
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
107 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
108 ## @example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
109 ## @group
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
110 ## n = 20;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
111 ## A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n)) + ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
112 ## toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
113 ## sparse (1, 2, 1, 1, n) * n / 2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
114 ## 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:
diff changeset
115 ## restart = 5;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
116 ## [M1, M2] = ilu (A); # in this tridiag case it corresponds to chol (A)'
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
117 ## M = M1 * M2;
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
118 ## Afcn = @@(x) A * x;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
119 ## Mfcn = @@(x) M \ x;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
120 ## M1fcn = @@(x) M1 \ x;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
121 ## M2fcn = @@(x) 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:
diff changeset
122 ## @end group
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
123 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
124 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
125 ## @sc{Example 1:} simplest usage of @code{tfqmr}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
126 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
127 ## @example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
128 ## x = tfqmr (A, b, [], n)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
129 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
130 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
131 ## @sc{Example 2:} @code{tfqmr} with a function which computes
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
132 ## @code{@var{A} * @var{x}}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
133 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
134 ## @example
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
135 ## x = tfqmr (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:
diff changeset
136 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
137 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
138 ## @sc{Example 3:} @code{tfqmr} 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:
diff changeset
139 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
140 ## @example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
141 ## x = tfqmr (A, b, [], 1e-06, n, M)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
142 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
143 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
144 ## @sc{Example 4:} @code{tfqmr} with a function as preconditioner
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
145 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
146 ## @example
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
147 ## x = tfqmr (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:
diff changeset
148 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
149 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
150 ## @sc{Example 5:} @code{tfqmr} 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:
diff changeset
151 ## 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:
diff changeset
152 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
153 ## @example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
154 ## x = tfqmr (A, b, [], 1e-6, n, M1, M2)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
155 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
156 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
157 ## @sc{Example 6:} @code{tfmqr} 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:
diff changeset
158 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
159 ## @example
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
160 ## x = tfqmr (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:
diff changeset
161 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
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:
diff changeset
163 ## @sc{Example 7:} @code{tfqmr} with as input a function requiring an argument
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
164 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
165 ## @example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
166 ## @group
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
167 ## function y = Ap (A, x, z) # compute A^z * x
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
168 ## y = x;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
169 ## for i = 1:z
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
170 ## y = A * y;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
171 ## endfor
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
172 ## endfunction
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
173 ## 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
174 ## x = tfqmr (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:
diff changeset
175 ## @end group
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
176 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
177 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
178 ## @sc{Example 8:} explicit example to show that @code{tfqmr} uses a
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
179 ## right preconditioner
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
180 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
181 ## @example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
182 ## @group
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
183 ## [M1, M2] = ilu (A + 0.3 * eye (n)); # factorization of A perturbed
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
184 ## M = M1 * M2;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
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:
diff changeset
186 ## ## reference solution computed by tfqmr after one iteration
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
187 ## [x_ref, fl] = tfqmr (A, b, [], 1, M)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
188 ##
25003
2365c2661b3c doc: Spellcheck documentation ahead of 4.4 release.
Rik <rik@octave.org>
parents: 24992
diff changeset
189 ## ## right preconditioning
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
190 ## [y, fl] = tfqmr (A / M, b, [], 1)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
191 ## x = M \ y # compare x and x_ref
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
192 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
193 ## @end group
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
194 ## @end example
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
195 ##
27984
b09432b20a84 maint: Remove special cases of old version control keywords in code base.
Rik <rik@octave.org>
parents: 27923
diff changeset
196 ## Reference:
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
197 ##
27984
b09432b20a84 maint: Remove special cases of old version control keywords in code base.
Rik <rik@octave.org>
parents: 27923
diff changeset
198 ## @nospell{Y. Saad}, @cite{Iterative Methods for Sparse Linear Systems},
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
199 ## Second edition, 2003, SIAM
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
200 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
201 ## @seealso{bicg, bicgstab, cgs, gmres, pcg, qmr, pcr}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
202 ##
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
203 ## @end deftypefn
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
204
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
205 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:
diff changeset
206 tfqmr (A, b, tol = [], maxit = [], M1 = [], M2 = [], ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
207 x0 = [], varargin)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
208
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
209 [Afcn, M1fcn, M2fcn] = __alltohandles__ (A, b, M1, M2, "tfqmr");
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
210
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
211 [tol, maxit, x0] = __default__input__ ({1e-06, 2 * min(20, rows (b)), ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
212 zeros(rows (b), 1)}, tol, ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
213 maxit, x0);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
214
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
215 maxit = 2 * maxit; # To be compatible, since iteration = odd+even ones
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
216
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
217 norm_b = norm (b, 2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
218 if (norm_b == 0)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
219 if (nargout < 2)
28912
0de38a6ef693 maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents: 27984
diff changeset
220 printf ("The right hand side vector is all zero so tfqmr \n")
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
221 printf ("returned an all zero solution without iterating.\n")
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
222 endif
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
223 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:
diff changeset
224 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:
diff changeset
225 flag = 0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
226 resvec = 0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
227 relres = 0;
28917
72d57dbcc305 maint: Add semicolon after break and return keywords.
Rik <rik@octave.org>
parents: 28912
diff changeset
228 return;
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
229 endif
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
230
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
231 x = x_pr = x_min = x0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
232 iter = iter_min = m = 0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
233 resvec = zeros (maxit, 1);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
234 flag = 1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
235
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
236 w = u = r = r_star = b - feval (Afcn, x0, varargin{:});
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
237 rho_1 = (r_star' * r);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
238 d = 0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
239 tau = norm (r, 2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
240 theta = eta = 0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
241 resvec (1, 1) = norm (r, 2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
242 it = 1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
243
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
244 try
28912
0de38a6ef693 maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents: 27984
diff changeset
245 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
246 u_hat = feval (M1fcn, u, varargin{:});
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
247 u_hat = feval (M2fcn, u_hat, varargin{:});
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
248 v = feval (Afcn, u_hat, varargin{:});
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
249 catch
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
250 flag = 2;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
251 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:
diff changeset
252 while ((flag != 2) && (iter < maxit) && ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
253 (resvec (iter + 1, 1) >= norm_b * tol))
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
254 if (it > 0) # iter is even
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
255 v_r = r_star' * v; # inner prod between r_star and v
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
256 if (v_r == 0)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
257 ## Essentially the next iteration doesn't change x,
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
258 ## and the iter after this will have a 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:
diff changeset
259 flag = 4;
28917
72d57dbcc305 maint: Add semicolon after break and return keywords.
Rik <rik@octave.org>
parents: 28912
diff changeset
260 break;
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
261 endif
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
262 alpha = rho_1 / v_r;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
263 u_1 = u - alpha * v; # u at the after iteration
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
264 endif
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
265 u_hat = feval (M1fcn, u, varargin{:});
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
266 u_hat = feval (M2fcn, u_hat, varargin{:});
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
267 w -= alpha * feval (Afcn, u_hat, varargin{:});
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
268 d = u_hat + ((theta * theta) / alpha) * eta * d;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
269 theta = norm (w, 2) / tau;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
270 c = 1 / sqrt (1 + theta * theta);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
271 tau *= theta * c;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
272 eta = (c * c) * alpha;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
273 x += eta * d;
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
274 r -= eta * feval (Afcn, d, varargin{:});
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
275 if (it < 0) # iter is odd
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
276 rho_2 = rho_1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
277 rho_1 = (r_star' * w);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
278 if (rho_1 == 0)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
279 ## Essentially the next iteration doesn't change x,
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
280 ## and the iter after this will have a 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:
diff changeset
281 flag = 4;
28917
72d57dbcc305 maint: Add semicolon after break and return keywords.
Rik <rik@octave.org>
parents: 28912
diff changeset
282 break;
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
283 endif
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
284 beta = rho_1 / rho_2;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
285 u_1 = w + beta * u; # u at the after iteration
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
286 u1_hat = feval (M1fcn, u_1, varargin{:});
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
287 u1_hat = feval (M2fcn, u1_hat, varargin{:});
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
288 v = feval (Afcn, u1_hat, varargin{:}) + ...
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
289 beta * (feval (Afcn, u_hat, varargin{:}) + beta * v);
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
290 endif
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
291 u = u_1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
292 iter += 1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
293 resvec (iter + 1, 1) = norm (r, 2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
294 if (resvec (iter + 1, 1) <= resvec (iter_min + 1, 1))
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
295 ## iter with min residual
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
296 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:
diff changeset
297 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:
diff changeset
298 endif
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
299 if (norm (x_pr - x) <= 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:
diff changeset
300 flag = 3; # Stagnation
28917
72d57dbcc305 maint: Add semicolon after break and return keywords.
Rik <rik@octave.org>
parents: 28912
diff changeset
301 break;
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
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:
diff changeset
303 x_pr = x;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
304 it = -it;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
305 endwhile
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
306 resvec = resvec (1: (iter + 1));
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
307
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
308 relres = resvec (iter_min + 1) / norm (b);
28912
0de38a6ef693 maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents: 27984
diff changeset
309 iter_min = floor (iter_min / 2); # compatibility, since it
0de38a6ef693 maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents: 27984
diff changeset
310 # makes two times the effective iterations
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
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:
diff changeset
312 if (relres <= tol)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
313 flag = 0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
314 endif
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
315
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
316 if (nargout < 2) # Output strings
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
317 switch (flag)
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
318 case {0}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
319 printf ("tfqmr 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:
diff changeset
320 printf ("to a solution with 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:
diff changeset
321 case {1}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
322 printf ("tfqmr 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:
diff changeset
323 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:
diff changeset
324 printf ("because the maximum number of iterations was reached.\n");
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
325 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:
diff changeset
326 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:
diff changeset
327 case {2}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
328 printf ("tfqmr 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:
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:
diff changeset
330 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:
diff changeset
331 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:
diff changeset
332 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:
diff changeset
333 case {3}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
334 printf ("tfqmr 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:
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:
diff changeset
336 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:
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:
diff changeset
338 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:
diff changeset
339 case {4}
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
340 printf ("tfqmr 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:
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:
diff changeset
342 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:
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:
diff changeset
344 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:
diff changeset
345 endswitch
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
346 endif
28945
6e460773bdda maint: Use newlines after "function" and before "endfunction" for clarity.
Rik <rik@octave.org>
parents: 28931
diff changeset
347
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
348 endfunction
28931
c5b1bbb95a66 maint: start %!demo or %!test blocks 2 newlines after endfunction.
Rik <rik@octave.org>
parents: 28929
diff changeset
349
c5b1bbb95a66 maint: start %!demo or %!test blocks 2 newlines after endfunction.
Rik <rik@octave.org>
parents: 28929
diff changeset
350
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
351 %!test
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
352 %! ## Check that all type of inputs work
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
353 %! A = toeplitz (sparse ([2, 1, 0, 0, 0]), sparse ([2, -1, 0, 0, 0]));
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
354 %! 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:
diff changeset
355 %! 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:
diff changeset
356 %! M2 = M1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
357 %! maxit = 10;
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
358 %! Afcn = @(z) A * z;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
359 %! M1_fcn = @(z) M1 \ z;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
360 %! M2_fcn = @(z) 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:
diff changeset
361 %! [x, flag] = tfqmr (A,b);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
362 %! 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:
diff changeset
363 %! [x, flag] = tfqmr (A, b, [], maxit, M1, M2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
364 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
365 %! [x, flag] = tfqmr (A, b, [], maxit, 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:
diff changeset
366 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
367 %! [x, flag] = tfqmr (A, b, [], maxit, 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:
diff changeset
368 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
369 %! [x, flag] = tfqmr (A, b, [], maxit, 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:
diff changeset
370 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
371 %! [x, flag] = tfqmr (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:
diff changeset
372 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
373 %! [x, flag] = tfqmr (Afcn, b, [], maxit, 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:
diff changeset
374 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
375 %! [x, flag] = tfqmr (Afcn, b, [], maxit, 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:
diff changeset
376 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
377 %! [x, flag] = tfqmr (Afcn, b, [], maxit, 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:
diff changeset
378 %! assert (flag, 0);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
379 %! [x, flag] = tfqmr (Afcn, b, [], maxit, 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:
diff changeset
380 %! 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:
diff changeset
381
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
382 %!shared A, b, n, M1, M2
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
383 %!
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
384 %!test
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
385 %! n = 100;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
386 %! A = spdiags ([-2*ones(n,1) 4*ones(n,1) -ones(n,1)], -1:1, n, n);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
387 %! 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:
diff changeset
388 %! tol = 1e-8;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
389 %! maxit = 15;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
390 %! M1 = spdiags ([ones(n,1)/(-2) ones(n,1)],-1:0, n, n);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
391 %! M2 = spdiags ([4*ones(n,1) -ones(n,1)], 0:1, n, n);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
392 %! [x, flag, relres, iter, resvec] = tfqmr (A, b, tol, maxit, M1, M2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
393 %! assert (x, ones (size (b)), 1e-7);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
394 %!
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
395 %!test
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
396 %!function y = afcn (x, a)
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
397 %! y = a * x;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
398 %!endfunction
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
399 %!
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
400 %! tol = 1e-8;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
401 %! maxit = 15;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
402 %!
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
403 %! [x, flag, relres, iter, resvec] = tfqmr (@(x) afcn (x, A), b,
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
404 %! tol, maxit, 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:
diff changeset
405 %! assert (x, ones (size (b)), 1e-7);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
406
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
407 %!test
24951
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
408 %! ## Jacobi preconditioner works
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
409 %! n = 10;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
410 %! tol = 1e-8;
24951
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
411 %! A = hilb (n) + 1i * hilb (n);
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
412 %! A(1,1) = 100;
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
413 %! A(n, 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:
diff changeset
414 %! b = sum (A, 2);
24951
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
415 %! [x, flag, relres, iter, resvec] = tfqmr (A, b, tol);
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
416 %! assert (x, ones (size (b)), 0.005);
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
417 %! assert (iter, 8);
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
418 %! [x, flag, relres, iter, resvec] = tfqmr (A, b, tol, [], diag (diag (A)));
24951
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
419 %! assert (x, ones (size (b)), 0.002);
27d68d7a482c Change a test in tfqmr to prevent random failures (bug #53319).
Marco Caliari <marco.caliari@univr.it>
parents: 24892
diff changeset
420 %! assert (iter, 6);
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
421
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
422 %!test
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
423 %! ## Solve complex linear system
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
424 %! A = [1 + 1i, 1 + 1i; 2 - 1i, 2 + 1i];
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
425 %! b = A * [1; 1];
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
426 %! [x, flag, relres, iter, resvec] = tfqmr (A, b, [], 3);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
427 %! assert (x, [1; 1], 1e-6);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
428
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
429 %!test
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
430 %! A = diag (1:50);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
431 %! A (1,50) = 10000;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
432 %! b = ones (50,1);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
433 %! [x, flag, relres, iter, resvec] = tfqmr (A, b, [], 100);
28929
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
434 %! assert (flag, 0);
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
435 %! assert (x, A \ b, 1e-05);
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
436 %! ## Detects a singular preconditioner
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
437 %! M = ones (50);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
438 %! M(1, 1) = 0;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
439 %! [x, flag] = tfqmr (A, b, [], 100, M);
28929
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
440 %! 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:
diff changeset
441
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
442 %!test
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
443 %! 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:
diff changeset
444 %! b = 1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
445 %! [x, flag] = tfqmr (A, b);
28929
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
446 %! 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:
diff changeset
447
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
448 %!test
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
449 %! A = 1;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
450 %! 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:
diff changeset
451 %! [x, flag] = tfqmr (A, b);
28929
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
452 %! 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:
diff changeset
453
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
454 %!test
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
455 %! 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:
diff changeset
456 %! 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:
diff changeset
457 %! [x, flag] = tfqmr (A, b);
28929
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
458 %! 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:
diff changeset
459
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
460 %!test
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
461 %!function y = Afcn (x)
28912
0de38a6ef693 maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents: 27984
diff changeset
462 %! A = toeplitz ([2, 1, 0, 0], [2, -1, 0, 0]);
0de38a6ef693 maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents: 27984
diff changeset
463 %! y = 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:
diff changeset
464 %!endfunction
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
465 %! [x, flag] = tfqmr ("Afcn", [1; 2; 2; 3]);
28929
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
466 %! 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:
diff changeset
467
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
468 %!test # unpreconditioned residual
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
469 %! A = toeplitz (sparse ([2, 1, 0, 0, 0]), sparse ([2, -1, 0, 0, 0]));
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
470 %! 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:
diff changeset
471 %! 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:
diff changeset
472 %! [x, flag, relres] = tfqmr (A, b, [], 3, M);
28929
9e43deb9bfc3 maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents: 28917
diff changeset
473 %! assert (relres, norm (b - A * x) / norm (b), 8 * eps);
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
474
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
475 %!demo # simplest use
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
476 %! n = 20;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
477 %! A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n)) + ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
478 %! toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
479 %! sparse (1, 2, 1, 1, n) * n / 2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
480 %! 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:
diff changeset
481 %! [M1, M2] = ilu (A + 0.1 * eye (n));
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
482 %! M = M1 * M2;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
483 %! x = tfqmr (A, b, [], n);
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
484 %! Afcn = @(x) A * x;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
485 %! x = tfqmr (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:
diff changeset
486 %! x = tfqmr (A, b, 1e-6, n, M);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
487 %! x = tfqmr (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
488 %! Mfcn = @(z) M \ z;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
489 %! x = tfqmr (Afcn, b, 1e-6, n, Mfcn);
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
490 %! M1fcn = @(z) M1 \ z;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
491 %! M2fcn = @(z) M2 \ z;
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
492 %! x = tfqmr (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:
diff changeset
493 %! function y = Ap (A, x, z) # compute A^z * x or (A^z)' * x
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
494 %! y = x;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
495 %! for i = 1:z
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
496 %! y = A * y;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
497 %! endfor
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
498 %! endfunction
30893
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
499 %! Afcn = @(x, p) Ap (A, x, p);
e1788b1a315f maint: Use "fcn" as preferred abbreviation for "function" in m-files.
Rik <rik@octave.org>
parents: 30564
diff changeset
500 %! x = tfqmr (Afcn, b, [], 2*n, [], [], [], 2); # solution of A^2 * x = b
24860
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
501
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
502 %!demo
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
503 %! n = 10;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
504 %! A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n)) + ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
505 %! toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
506 %! sparse (1, 2, 1, 1, n) * n / 2);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
507 %! 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:
diff changeset
508 %! [M1, M2] = ilu (A + 0.3 * eye (n)); # factorization of A perturbed
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
509 %! M = M1 * M2;
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
510 %!
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
511 %! ## reference solution computed by tfqmr after one iteration
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
512 %! [x_ref, fl] = tfqmr (A, b, [], 1, M);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
513 %! x_ref
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
514 %!
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
515 %! ## right preconditioning
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
516 %! [y, fl] = tfqmr (A / M, b, [], 1);
6266e321ef22 Merge code from SoCiS 2016 "Improve iterative methods" project (http://socis16octave-improveiterativemethods.blogspot.com)
Cristiano Dorigo <cristiano.dorigo@hotmail.it>
parents:
diff changeset
517 %! x = M \ y # compare x and x_ref