annotate libinterp/corefcn/qz.cc @ 33584:3fe954c2fd25 default tip @

maint: merge stable to default
author Rik <rik@octave.org>
date Mon, 13 May 2024 11:41:11 -0700
parents f53ac65ffba6
children
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 //
32632
2e484f9f1f18 maint: update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents: 32488
diff changeset
3 // Copyright (C) 1998-2024 The Octave Project Developers
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
4 //
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/>.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
7 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
8 // This file is part of Octave.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
9 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
10 // Octave is free software: you can redistribute it and/or modify it
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
11 // under the terms of the GNU General Public License as published by
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
12 // the Free Software Foundation, either version 3 of the License, or
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
13 // (at your option) any later version.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
14 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
15 // Octave is distributed in the hope that it will be useful, but
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
18 // GNU General Public License for more details.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
19 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
20 // You should have received a copy of the GNU General Public License
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
21 // along with Octave; see the file COPYING. If not, see
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
22 // <https://www.gnu.org/licenses/>.
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 ////////////////////////////////////////////////////////////////////////
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
25
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
26 // Generalized eigenvalue balancing via LAPACK
3911
8389e78e67d4 [project @ 2002-04-28 02:15:38 by jwe]
jwe
parents: 3887
diff changeset
27
25570
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
28 // Originally written by A. S. Hodel <scotte@eng.auburn.edu>, but is
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
29 // substantially different with the change to use LAPACK.
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
30
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
31 #undef DEBUG
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
32 #undef DEBUG_SORT
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
33 #undef DEBUG_EIG
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
34
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
35 #if defined (HAVE_CONFIG_H)
21301
40de9f8f23a6 Use '#include "config.h"' rather than <config.h>.
Rik <rik@octave.org>
parents: 21279
diff changeset
36 # include "config.h"
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 9758
diff changeset
37 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
38
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
39 #include <cctype>
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
40 #include <cmath>
4051
b79da8779a0e [project @ 2002-08-17 19:38:32 by jwe]
jwe
parents: 3911
diff changeset
41
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
42 #if defined (DEBUG_EIG)
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
43 # include <iomanip>
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
44 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
45
4153
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 4051
diff changeset
46 #include "f77-fcn.h"
22322
93b3cdd36854 move most f77 function decls to separate header files
John W. Eaton <jwe@octave.org>
parents: 22317
diff changeset
47 #include "lo-lapack-proto.h"
21279
eb1524b07fe3 better use of templates for qr classes
John W. Eaton <jwe@octave.org>
parents: 21200
diff changeset
48 #include "qr.h"
4153
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 4051
diff changeset
49 #include "quit.h"
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 4051
diff changeset
50
15039
e753177cde93 maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents: 14876
diff changeset
51 #include "defun.h"
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
52 #include "error.h"
21100
e39e05d90788 Switch gripe_XXX to either err_XXX or warn_XXX naming scheme.
Rik <rik@octave.org>
parents: 21024
diff changeset
53 #include "errwarn.h"
20940
48b2ad5ee801 maint: Rename oct-obj.[cc|h] to ovl.[cc|h] for clarity.
Rik <rik@octave.org>
parents: 20892
diff changeset
54 #include "ovl.h"
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
55 #if defined (DEBUG) || defined (DEBUG_SORT)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
56 # include "pager.h"
21200
fcac5dbbf9ed maint: Indent #ifdef blocks in libinterp.
Rik <rik@octave.org>
parents: 21110
diff changeset
57 # include "pr-output.h"
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
58 #endif
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
59
31605
e88a07dec498 maint: Use macros to begin/end C++ namespaces.
Rik <rik@octave.org>
parents: 30564
diff changeset
60 OCTAVE_BEGIN_NAMESPACE(octave)
29958
32c3a5805893 move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents: 29558
diff changeset
61
15039
e753177cde93 maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents: 14876
diff changeset
62 DEFUN (qz, args, nargout,
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
63 doc: /* -*- texinfo -*-
32488
d7a3ed7f2fdc doc: grammarcheck C++ files before 9.1 release.
Rik <rik@octave.org>
parents: 31977
diff changeset
64 @deftypefn {} {[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}] =} qz (@var{A}, @var{B})
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
65 @deftypefnx {} {[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}] =} qz (@var{A}, @var{B}, @var{opt})
25106
d7ad543255c5 doc: Shorten very long first sentences of docstrings (bug #53388).
Rik <rik@octave.org>
parents: 25103
diff changeset
66 Compute the QZ@tie{}decomposition of a generalized eigenvalue problem.
d7ad543255c5 doc: Shorten very long first sentences of docstrings (bug #53388).
Rik <rik@octave.org>
parents: 25103
diff changeset
67
d7ad543255c5 doc: Shorten very long first sentences of docstrings (bug #53388).
Rik <rik@octave.org>
parents: 25103
diff changeset
68 The generalized eigenvalue problem is defined as
d7ad543255c5 doc: Shorten very long first sentences of docstrings (bug #53388).
Rik <rik@octave.org>
parents: 25103
diff changeset
69
23428
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
70 @tex
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
71 $$A x = \lambda B x$$
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
72 @end tex
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
73 @ifnottex
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
74
23428
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
75 @math{A x = @var{lambda} B x}
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
76
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
77 @end ifnottex
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
78
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
79 There are two calling forms of the function:
23428
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
80
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
81 @enumerate
23428
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
82 @item @code{[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}, @var{lambda}] = qz (@var{A}, @var{B})}
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
83
32488
d7a3ed7f2fdc doc: grammarcheck C++ files before 9.1 release.
Rik <rik@octave.org>
parents: 31977
diff changeset
84 Compute the complex QZ@tie{}decomposition, generalized eigenvectors, and
d7a3ed7f2fdc doc: grammarcheck C++ files before 9.1 release.
Rik <rik@octave.org>
parents: 31977
diff changeset
85 generalized eigenvalues.
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
86 @tex
29556
d75aa2bf4915 qz.cc: return correct number of eigenvalues (bug #60357).
Rik <rik@octave.org>
parents: 29358
diff changeset
87 $$ AA = Q^T AZ, BB = Q^T BZ $$
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
88 $$ { \rm diag }(BB)AV = BV{ \rm diag }(AA) $$
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
89 $$ { \rm diag }(BB) W^T A = { \rm diag }(AA) W^T B $$
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
90 @end tex
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
91 @ifnottex
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
92
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
93 @example
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
94 @group
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
95
29556
d75aa2bf4915 qz.cc: return correct number of eigenvalues (bug #60357).
Rik <rik@octave.org>
parents: 29358
diff changeset
96 @var{AA} = @var{Q} * @var{A} * @var{Z}, @var{BB} = @var{Q} * @var{B} * @var{Z}
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
97 @var{A} * @var{V} * diag (diag (@var{BB})) = @var{B} * @var{V} * diag (diag (@var{AA}))
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
98 diag (diag (@var{BB})) * @var{W}' * @var{A} = diag (diag (@var{AA})) * @var{W}' * @var{B}
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
99
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
100 @end group
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
101 @end example
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
102
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
103 @end ifnottex
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
104 with @var{AA} and @var{BB} upper triangular, and @var{Q} and @var{Z}
32488
d7a3ed7f2fdc doc: grammarcheck C++ files before 9.1 release.
Rik <rik@octave.org>
parents: 31977
diff changeset
105 unitary. The matrices @var{V} and @var{W} respectively contain the right
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
106 and left generalized eigenvectors.
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
107
23428
1bc0e610e293 doc: Redo docstring for qz (bug #50846).
Rik <rik@octave.org>
parents: 23219
diff changeset
108 @item @code{[@var{AA}, @var{BB}, @var{Z} @{, @var{lambda}@}] = qz (@var{A}, @var{B}, @var{opt})}
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
109
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
110 The @var{opt} argument must be equal to either @qcode{"real"} or
32488
d7a3ed7f2fdc doc: grammarcheck C++ files before 9.1 release.
Rik <rik@octave.org>
parents: 31977
diff changeset
111 @qcode{"complex"}. If it is equal to @qcode{"complex"}, then this
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
112 calling form is equivalent to the first one with only two input
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
113 arguments.
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
114
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
115 If @var{opt} is equal to @qcode{"real"}, then the real QZ decomposition
32488
d7a3ed7f2fdc doc: grammarcheck C++ files before 9.1 release.
Rik <rik@octave.org>
parents: 31977
diff changeset
116 is computed. In particular, @var{AA} is only guaranteed to be
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
117 quasi-upper triangular with 1-by-1 and 2-by-2 blocks on the diagonal,
32488
d7a3ed7f2fdc doc: grammarcheck C++ files before 9.1 release.
Rik <rik@octave.org>
parents: 31977
diff changeset
118 and @var{Q} and @var{Z} are orthogonal. The identities mentioned above
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
119 for right and left generalized eigenvectors are only verified if
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
120 @var{AA} is upper triangular (i.e., when all the generalized eigenvalues
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
121 are real, in which case the real and complex QZ coincide).
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
122
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
123 @end enumerate
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
124
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
125 Note: @code{qz} performs permutation balancing, but not scaling
28961
d9d028b479ac doc: Use @code{} within alternate text for @xref,@pxref macros in libinterp/
Rik <rik@octave.org>
parents: 27923
diff changeset
126 (@pxref{XREFbalance,,@code{balance}}), which may be lead to less accurate
d9d028b479ac doc: Use @code{} within alternate text for @xref,@pxref macros in libinterp/
Rik <rik@octave.org>
parents: 27923
diff changeset
127 results than @code{eig}. The order of output arguments was selected for
d9d028b479ac doc: Use @code{} within alternate text for @xref,@pxref macros in libinterp/
Rik <rik@octave.org>
parents: 27923
diff changeset
128 compatibility with @sc{matlab}.
29556
d75aa2bf4915 qz.cc: return correct number of eigenvalues (bug #60357).
Rik <rik@octave.org>
parents: 29358
diff changeset
129 @seealso{eig, gsvd, balance, chol, hess, lu, qr, qzhess, schur}
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21751
diff changeset
130 @end deftypefn */)
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
131 {
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
132 int nargin = args.length ();
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
133
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
134 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
135 octave_stdout << "qz: nargin = " << nargin
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
136 << ", nargout = " << nargout << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
137 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
138
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
139 if (nargin < 2 || nargin > 3 || nargout > 6)
20802
8bb38ba1bad6 eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents: 20700
diff changeset
140 print_usage ();
8bb38ba1bad6 eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents: 20700
diff changeset
141
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
142 bool complex_case = true;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
143
24116
e1020353fff5 qz: initialize optional ordering flag safely
Mike Miller <mtmiller@octave.org>
parents: 23826
diff changeset
144 if (nargin == 3)
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
145 {
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
146 std::string opt = args(2).xstring_value ("qz: OPT must be a string");
8927
d75f4ee0538d qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
147
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
148 if (opt == "real")
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
149 {
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
150 if (args(0).iscomplex () || args(1).iscomplex ())
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
151 {
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
152 warning ("qz: ignoring 'real' option with complex matrices");
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
153 complex_case = true;
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
154 }
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
155 else
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
156 complex_case = false;
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
157 }
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
158 else if (opt == "complex")
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
159 complex_case = true;
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
160 else
31977
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
161 error ("qz: OPT must be 'real' or 'complex'");
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
162 }
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
163
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
164 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
165 octave_stdout << "qz: check matrix A" << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
166 #endif
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
167
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
168 // Matrix A: check dimensions.
29961
7d6709900da7 eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents: 29958
diff changeset
169 F77_INT nn = to_f77_int (args(0).rows ());
7d6709900da7 eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents: 29958
diff changeset
170 F77_INT nc = to_f77_int (args(0).columns ());
20892
c07bee629973 2015 Code Sprint: use ovl ().
Rik <rik@octave.org>
parents: 20853
diff changeset
171
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
172 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
173 octave_stdout << "Matrix A dimensions: (" << nn << ',' << nc << ')'
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
174 << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
175 #endif
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
176
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
177 if (nc != nn)
21110
3d0d84305600 Use err_square_matrix_required more widely.
Rik <rik@octave.org>
parents: 21100
diff changeset
178 err_square_matrix_required ("qz", "A");
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
179
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
180 // Matrix A: get value.
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
181 Matrix aa;
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
182 ComplexMatrix caa;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
183
23581
c3075ae020e1 maint: Deprecate is_complex_type and replace with iscomplex.
Rik <rik@octave.org>
parents: 23577
diff changeset
184 if (args(0).iscomplex ())
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
185 caa = args(0).complex_matrix_value ();
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
186 else
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
187 aa = args(0).matrix_value ();
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
188
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
189 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
190 octave_stdout << "qz: check matrix B" << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
191 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
192
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
193 // Extract argument 2 (bb, or cbb if complex).
29961
7d6709900da7 eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents: 29958
diff changeset
194 F77_INT b_nr = to_f77_int (args(1).rows ());
7d6709900da7 eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents: 29958
diff changeset
195 F77_INT b_nc = to_f77_int (args(1).columns ());
22964
0c12642be005 use F77_INT instead of octave_idx_type for libinterp QZ function
John W. Eaton <jwe@octave.org>
parents: 22860
diff changeset
196
0c12642be005 use F77_INT instead of octave_idx_type for libinterp QZ function
John W. Eaton <jwe@octave.org>
parents: 22860
diff changeset
197 if (nn != b_nc || nn != b_nr)
29958
32c3a5805893 move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents: 29558
diff changeset
198 ::err_nonconformant ();
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
199
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
200 Matrix bb;
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
201 ComplexMatrix cbb;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
202
23581
c3075ae020e1 maint: Deprecate is_complex_type and replace with iscomplex.
Rik <rik@octave.org>
parents: 23577
diff changeset
203 if (args(1).iscomplex ())
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
204 cbb = args(1).complex_matrix_value ();
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
205 else
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
206 bb = args(1).matrix_value ();
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
207
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
208 // First, declare variables used in both the real and complex cases.
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
209 // FIXME: There are a lot of excess variables declared.
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
210 // Probably a better way to handle this.
30390
a61e1a0f6024 maint: style check C++ files in libinterp/ ahead of 7.1 release.
Rik <rik@octave.org>
parents: 29961
diff changeset
211 Matrix QQ (nn, nn), ZZ (nn, nn), VR (nn, nn), VL (nn, nn);
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
212 RowVector alphar (nn), alphai (nn), betar (nn);
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
213 ComplexRowVector xalpha (nn), xbeta (nn);
30390
a61e1a0f6024 maint: style check C++ files in libinterp/ ahead of 7.1 release.
Rik <rik@octave.org>
parents: 29961
diff changeset
214 ComplexMatrix CQ (nn, nn), CZ (nn, nn), CVR (nn, nn), CVL (nn, nn);
22964
0c12642be005 use F77_INT instead of octave_idx_type for libinterp QZ function
John W. Eaton <jwe@octave.org>
parents: 22860
diff changeset
215 F77_INT ilo, ihi, info;
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
216 char comp_q = (nargout >= 3 ? 'V' : 'N');
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
217 char comp_z = (nargout >= 4 ? 'V' : 'N');
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
218
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
219 // Initialize Q, Z to identity matrix if either is needed
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
220 if (comp_q == 'V' || comp_z == 'V')
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
221 {
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
222 double *QQptr = QQ.rwdata ();
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
223 double *ZZptr = ZZ.rwdata ();
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
224 std::fill_n (QQptr, QQ.numel (), 0.0);
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
225 std::fill_n (ZZptr, ZZ.numel (), 0.0);
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
226 for (F77_INT i = 0; i < nn; i++)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
227 {
30390
a61e1a0f6024 maint: style check C++ files in libinterp/ ahead of 7.1 release.
Rik <rik@octave.org>
parents: 29961
diff changeset
228 QQ(i, i) = 1.0;
a61e1a0f6024 maint: style check C++ files in libinterp/ ahead of 7.1 release.
Rik <rik@octave.org>
parents: 29961
diff changeset
229 ZZ(i, i) = 1.0;
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
230 }
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
231 }
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
232
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
233 // Always perform permutation balancing.
4552
6f3382e08a52 [project @ 2003-10-27 20:38:02 by jwe]
jwe
parents: 4153
diff changeset
234 const char bal_job = 'P';
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
235 RowVector lscale (nn), rscale (nn), work (6 * nn), rwork (nn);
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
236
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
237 if (complex_case)
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
238 {
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
239 #if defined (DEBUG)
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
240 if (comp_q == 'V')
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
241 octave_stdout << "qz: performing balancing; CQ =\n" << CQ << std::endl;
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
242 #endif
23582
0cc2011d800e maint: Deprecate is_real_type and replace with isreal.
Rik <rik@octave.org>
parents: 23581
diff changeset
243 if (args(0).isreal ())
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
244 caa = ComplexMatrix (aa);
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
245
23582
0cc2011d800e maint: Deprecate is_real_type and replace with isreal.
Rik <rik@octave.org>
parents: 23581
diff changeset
246 if (args(1).isreal ())
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
247 cbb = ComplexMatrix (bb);
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
248
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
249 if (comp_q == 'V')
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
250 CQ = ComplexMatrix (QQ);
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
251
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
252 if (comp_z == 'V')
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
253 CZ = ComplexMatrix (ZZ);
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
254
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
255 F77_XFCN (zggbal, ZGGBAL,
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
256 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
257 nn, F77_DBLE_CMPLX_ARG (caa.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
258 F77_DBLE_CMPLX_ARG (cbb.rwdata ()),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
259 nn, ilo, ihi, lscale.rwdata (),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
260 rscale.rwdata (), work.rwdata (), info
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
261 F77_CHAR_ARG_LEN (1)));
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
262 }
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
263 else
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
264 {
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
265 #if defined (DEBUG)
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
266 if (comp_q == 'V')
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
267 octave_stdout << "qz: performing balancing; QQ =\n" << QQ << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
268 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
269
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
270 F77_XFCN (dggbal, DGGBAL,
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
271 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
272 nn, aa.rwdata (), nn, bb.rwdata (),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
273 nn, ilo, ihi, lscale.rwdata (),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
274 rscale.rwdata (), work.rwdata (), info
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
275 F77_CHAR_ARG_LEN (1)));
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
276 }
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
277
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
278 // Only permutation balance above is done. Skip scaling balance.
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
279
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
280 char qz_job = (nargout < 2 ? 'E' : 'S');
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
281
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
282 if (complex_case)
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
283 {
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
284 // Complex case.
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
285
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
286 // The QR decomposition of cbb.
29961
7d6709900da7 eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents: 29958
diff changeset
287 math::qr<ComplexMatrix> cbqr (cbb);
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
288 // The R matrix of QR decomposition for cbb.
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
289 cbb = cbqr.R ();
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
290 // (Q*)caa for following work.
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
291 caa = (cbqr.Q ().hermitian ()) * caa;
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
292 CQ = CQ * cbqr.Q ();
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
293
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
294 F77_XFCN (zgghrd, ZGGHRD,
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
295 (F77_CONST_CHAR_ARG2 (&comp_q, 1),
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
296 F77_CONST_CHAR_ARG2 (&comp_z, 1),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
297 nn, ilo, ihi, F77_DBLE_CMPLX_ARG (caa.rwdata ()),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
298 nn, F77_DBLE_CMPLX_ARG (cbb.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
299 F77_DBLE_CMPLX_ARG (CQ.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
300 F77_DBLE_CMPLX_ARG (CZ.rwdata ()), nn, info
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
301 F77_CHAR_ARG_LEN (1)
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
302 F77_CHAR_ARG_LEN (1)));
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
303
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
304 ComplexRowVector cwork (nn);
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
305
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
306 F77_XFCN (zhgeqz, ZHGEQZ,
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
307 (F77_CONST_CHAR_ARG2 (&qz_job, 1),
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
308 F77_CONST_CHAR_ARG2 (&comp_q, 1),
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
309 F77_CONST_CHAR_ARG2 (&comp_z, 1),
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
310 nn, ilo, ihi,
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
311 F77_DBLE_CMPLX_ARG (caa.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
312 F77_DBLE_CMPLX_ARG (cbb.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
313 F77_DBLE_CMPLX_ARG (xalpha.rwdata ()),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
314 F77_DBLE_CMPLX_ARG (xbeta.rwdata ()),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
315 F77_DBLE_CMPLX_ARG (CQ.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
316 F77_DBLE_CMPLX_ARG (CZ.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
317 F77_DBLE_CMPLX_ARG (cwork.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
318 rwork.rwdata (), info
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
319 F77_CHAR_ARG_LEN (1)
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
320 F77_CHAR_ARG_LEN (1)
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
321 F77_CHAR_ARG_LEN (1)));
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
322
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
323 if (comp_q == 'V')
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
324 {
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
325 // Left eigenvector.
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
326 F77_XFCN (zggbak, ZGGBAK,
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
327 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
328 F77_CONST_CHAR_ARG2 ("L", 1),
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
329 nn, ilo, ihi, lscale.data (), rscale.data (),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
330 nn, F77_DBLE_CMPLX_ARG (CQ.rwdata ()), nn, info
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
331 F77_CHAR_ARG_LEN (1)
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
332 F77_CHAR_ARG_LEN (1)));
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
333 }
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
334
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
335 if (comp_z == 'V')
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
336 {
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
337 // Right eigenvector.
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
338 F77_XFCN (zggbak, ZGGBAK,
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
339 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
340 F77_CONST_CHAR_ARG2 ("R", 1),
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
341 nn, ilo, ihi, lscale.data (), rscale.data (),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
342 nn, F77_DBLE_CMPLX_ARG (CZ.rwdata ()), nn, info
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
343 F77_CHAR_ARG_LEN (1)
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
344 F77_CHAR_ARG_LEN (1)));
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
345 }
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
346
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
347 }
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
348 else
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
349 {
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
350 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
351 octave_stdout << "qz: performing qr decomposition of bb" << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
352 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
353
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
354 // Compute the QR factorization of bb.
29961
7d6709900da7 eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents: 29958
diff changeset
355 math::qr<Matrix> bqr (bb);
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
356
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
357 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
358 octave_stdout << "qz: qr (bb) done; now performing qz decomposition"
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
359 << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
360 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
361
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
362 bb = bqr.R ();
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
363
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
364 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
365 octave_stdout << "qz: extracted bb" << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
366 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
367
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
368 aa = (bqr.Q ()).transpose () * aa;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
369
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
370 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
371 octave_stdout << "qz: updated aa " << std::endl;
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
372 octave_stdout << "bqr.Q () =\n" << bqr.Q () << std::endl;
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
373
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
374 if (comp_q == 'V')
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
375 octave_stdout << "QQ =" << QQ << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
376 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
377
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
378 if (comp_q == 'V')
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
379 QQ = QQ * bqr.Q ();
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
380
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
381 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
382 octave_stdout << "qz: precursors done..." << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
383 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
384
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
385 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
386 octave_stdout << "qz: comp_q = " << comp_q << ", comp_z = " << comp_z
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
387 << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
388 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
389
21562
6c2fd62db1f7 maint: Eliminate accidental double spaces in code.
Rik <rik@octave.org>
parents: 21301
diff changeset
390 // Reduce to generalized Hessenberg form.
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
391 F77_XFCN (dgghrd, DGGHRD,
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
392 (F77_CONST_CHAR_ARG2 (&comp_q, 1),
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
393 F77_CONST_CHAR_ARG2 (&comp_z, 1),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
394 nn, ilo, ihi, aa.rwdata (),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
395 nn, bb.rwdata (), nn, QQ.rwdata (), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
396 ZZ.rwdata (), nn, info
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
397 F77_CHAR_ARG_LEN (1)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
398 F77_CHAR_ARG_LEN (1)));
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
399
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
400 // Check if just computing generalized eigenvalues,
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
401 // or if we're actually computing the decomposition.
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
402
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
403 // Reduce to generalized Schur form.
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
404 F77_XFCN (dhgeqz, DHGEQZ,
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
405 (F77_CONST_CHAR_ARG2 (&qz_job, 1),
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
406 F77_CONST_CHAR_ARG2 (&comp_q, 1),
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
407 F77_CONST_CHAR_ARG2 (&comp_z, 1),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
408 nn, ilo, ihi, aa.rwdata (), nn, bb.rwdata (),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
409 nn, alphar.rwdata (), alphai.rwdata (),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
410 betar.rwdata (), QQ.rwdata (), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
411 ZZ.rwdata (), nn, work.rwdata (), nn, info
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
412 F77_CHAR_ARG_LEN (1)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
413 F77_CHAR_ARG_LEN (1)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
414 F77_CHAR_ARG_LEN (1)));
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
415
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
416 if (comp_q == 'V')
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
417 {
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
418 F77_XFCN (dggbak, DGGBAK,
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
419 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
420 F77_CONST_CHAR_ARG2 ("L", 1),
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
421 nn, ilo, ihi, lscale.data (), rscale.data (),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
422 nn, QQ.rwdata (), nn, info
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
423 F77_CHAR_ARG_LEN (1)
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
424 F77_CHAR_ARG_LEN (1)));
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
425
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
426 #if defined (DEBUG)
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
427 if (comp_q == 'V')
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
428 octave_stdout << "qz: balancing done; QQ=\n" << QQ << std::endl;
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
429 #endif
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
430 }
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
431
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
432 // then right
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
433 if (comp_z == 'V')
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
434 {
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
435 F77_XFCN (dggbak, DGGBAK,
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
436 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
437 F77_CONST_CHAR_ARG2 ("R", 1),
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
438 nn, ilo, ihi, lscale.data (), rscale.data (),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
439 nn, ZZ.rwdata (), nn, info
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
440 F77_CHAR_ARG_LEN (1)
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
441 F77_CHAR_ARG_LEN (1)));
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
442
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
443 #if defined (DEBUG)
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
444 if (comp_z == 'V')
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
445 octave_stdout << "qz: balancing done; ZZ=\n" << ZZ << std::endl;
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
446 #endif
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
447 }
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
448
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
449 }
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
450
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
451 // Right, left eigenvector matrices.
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
452 if (nargout >= 5)
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
453 {
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
454 // Which side to compute?
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
455 char side = (nargout == 5 ? 'R' : 'B');
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
456 // Compute all of them and backtransform
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
457 char howmany = 'B';
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
458 // Dummy pointer; select is not used.
23457
21baad6b35c4 maint: Use C++11 nullptr rather than 0 or NULL when possible.
Rik <rik@octave.org>
parents: 23449
diff changeset
459 F77_INT *select = nullptr;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
460
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
461 if (complex_case)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
462 {
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
463 CVL = CQ;
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
464 CVR = CZ;
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
465 ComplexRowVector cwork2 (2 * nn);
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
466 RowVector rwork2 (8 * nn);
22964
0c12642be005 use F77_INT instead of octave_idx_type for libinterp QZ function
John W. Eaton <jwe@octave.org>
parents: 22860
diff changeset
467 F77_INT m;
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
468
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
469 F77_XFCN (ztgevc, ZTGEVC,
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
470 (F77_CONST_CHAR_ARG2 (&side, 1),
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
471 F77_CONST_CHAR_ARG2 (&howmany, 1),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
472 select, nn, F77_DBLE_CMPLX_ARG (caa.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
473 F77_DBLE_CMPLX_ARG (cbb.rwdata ()),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
474 nn, F77_DBLE_CMPLX_ARG (CVL.rwdata ()), nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
475 F77_DBLE_CMPLX_ARG (CVR.rwdata ()), nn, nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
476 m, F77_DBLE_CMPLX_ARG (cwork2.rwdata ()),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
477 rwork2.rwdata (), info
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
478 F77_CHAR_ARG_LEN (1)
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
479 F77_CHAR_ARG_LEN (1)));
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
480 }
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
481 else
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
482 {
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
483 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
484 octave_stdout << "qz: computing generalized eigenvectors" << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
485 #endif
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
486
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
487 VL = QQ;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
488 VR = ZZ;
22964
0c12642be005 use F77_INT instead of octave_idx_type for libinterp QZ function
John W. Eaton <jwe@octave.org>
parents: 22860
diff changeset
489 F77_INT m;
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
490
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
491 F77_XFCN (dtgevc, DTGEVC,
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
492 (F77_CONST_CHAR_ARG2 (&side, 1),
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
493 F77_CONST_CHAR_ARG2 (&howmany, 1),
32660
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
494 select, nn, aa.rwdata (), nn, bb.rwdata (),
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
495 nn, VL.rwdata (), nn, VR.rwdata (), nn, nn,
f53ac65ffba6 maint: New method rwdata() as clearer alternative to fortran_vec().
Rik <rik@octave.org>
parents: 32632
diff changeset
496 m, work.rwdata (), info
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
497 F77_CHAR_ARG_LEN (1)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
498 F77_CHAR_ARG_LEN (1)));
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
499 }
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
500 }
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
501
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
502 octave_value_list retval (nargout);
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
503
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
504 switch (nargout)
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
505 {
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
506 case 6:
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
507 // Return left eigenvectors.
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
508 if (complex_case)
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
509 retval(5) = CVL;
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
510 else
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
511 retval(5) = VL;
23826
d69021d58a61 avoid fallthrough warnings
John W. Eaton <jwe@octave.org>
parents: 23807
diff changeset
512 OCTAVE_FALLTHROUGH;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
513
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
514 case 5:
23430
32d532b8e7d0 qz.cc: Overhaul qz function.
Rik <rik@octave.org>
parents: 23429
diff changeset
515 // Return right eigenvectors.
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
516 if (complex_case)
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
517 retval(4) = CVR;
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
518 else
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
519 retval(4) = VR;
23826
d69021d58a61 avoid fallthrough warnings
John W. Eaton <jwe@octave.org>
parents: 23807
diff changeset
520 OCTAVE_FALLTHROUGH;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
521
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
522 case 4:
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
523 if (complex_case)
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
524 retval(3) = CZ;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
525 else
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
526 retval(3) = ZZ;
23826
d69021d58a61 avoid fallthrough warnings
John W. Eaton <jwe@octave.org>
parents: 23807
diff changeset
527 OCTAVE_FALLTHROUGH;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
528
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
529 case 3:
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
530 if (complex_case)
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
531 retval(2) = CQ.hermitian ();
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
532 else
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
533 retval(2) = QQ.transpose ();
23826
d69021d58a61 avoid fallthrough warnings
John W. Eaton <jwe@octave.org>
parents: 23807
diff changeset
534 OCTAVE_FALLTHROUGH;
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
535
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
536 case 2:
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
537 case 1:
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
538 case 0:
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
539 {
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
540 if (complex_case)
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
541 {
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
542 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
543 octave_stdout << "qz: retval(1) = cbb =\n";
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
544 octave_print_internal (octave_stdout, cbb);
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
545 octave_stdout << "\nqz: retval(0) = caa =\n";
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
546 octave_print_internal (octave_stdout, caa);
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
547 octave_stdout << std::endl;
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
548 #endif
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
549 retval(1) = cbb;
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
550 retval(0) = caa;
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
551 }
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
552 else
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
553 {
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
554 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
555 octave_stdout << "qz: retval(1) = bb =\n";
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
556 octave_print_internal (octave_stdout, bb);
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
557 octave_stdout << "\nqz: retval(0) = aa =\n";
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
558 octave_print_internal (octave_stdout, aa);
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
559 octave_stdout << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
560 #endif
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
561 retval(1) = bb;
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
562 retval(0) = aa;
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
563 }
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
564 }
10311
a217e1d74353 untabify qz.cc
John W. Eaton <jwe@octave.org>
parents: 10307
diff changeset
565 break;
10307
4e4270ab70d6 support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents: 10155
diff changeset
566
25103
078b795c5219 maint: style check C++ ahead of 4.4 release.
Rik <rik@octave.org>
parents: 25054
diff changeset
567 }
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
568
31977
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
569 // FIXME: The API for qz changed in version 9.
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
570 // These warnings can be removed in Octave version 11.
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
571 if (nargout == 1)
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
572 {
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
573 warning_with_id ("Octave:qz:single-arg-out",
31977
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
574 "qz: requesting a single output argument no longer returns eigenvalues since version 9");
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
575 disable_warning ("Octave:qz:single-arg-out");
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
576 }
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
577
31977
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
578 if (nargin == 2 && args(0).isreal () && args(1).isreal ()
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
579 && retval(0).iscomplex ())
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
580 {
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
581 warning_with_id ("Octave:qz:complex-default",
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
582 "qz: returns the complex QZ by default on real matrices since version 9");
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
583 disable_warning ("Octave:qz:complex-default");
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
584 }
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
585
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21562
diff changeset
586 #if defined (DEBUG)
26001
acb4689aa5f2 * qz.cc: Use octave_stdout instead of std::cout for debugging messages.
John W. Eaton <jwe@octave.org>
parents: 25734
diff changeset
587 octave_stdout << "qz: exiting (at long last)" << std::endl;
3185
9580887dd160 [project @ 1998-09-26 02:45:55 by jwe]
jwe
parents: 3183
diff changeset
588 #endif
3183
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
589
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
590 return retval;
5edc539c2f80 [project @ 1998-09-25 03:37:22 by jwe]
jwe
parents:
diff changeset
591 }
13074
0b789a03bde1 codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12642
diff changeset
592
0b789a03bde1 codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12642
diff changeset
593 /*
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
594 ## Example 7.7.3 in Golub & Van Loan (3rd ed.; not present in 4th ed.)
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
595 ## The generalized eigenvalues are real, so the real and complex QZ coincide.
13074
0b789a03bde1 codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12642
diff changeset
596 %!test
0b789a03bde1 codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12642
diff changeset
597 %! a = [ 10 1 2;
0b789a03bde1 codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12642
diff changeset
598 %! 1 2 -1;
29556
d75aa2bf4915 qz.cc: return correct number of eigenvalues (bug #60357).
Rik <rik@octave.org>
parents: 29358
diff changeset
599 %! 1 1 2 ];
d75aa2bf4915 qz.cc: return correct number of eigenvalues (bug #60357).
Rik <rik@octave.org>
parents: 29358
diff changeset
600 %! b = reshape (1:9, 3,3);
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
601 %! [aa, bb, q, z, v, w] = qz (a, b);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
602 %! assert (isreal (aa) && isreal (bb) && isreal (q) && isreal (z) ...
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
603 %! && isreal (v) && isreal (w));
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
604 %! assert (istriu (aa) && istriu (bb));
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
605 %! assert (q * q', eye(3), 1e-15);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
606 %! assert (z * z', eye(3), 1e-15);
13088
4ffea87ad71b codesprint: Fix tolerance for qz.cc tests.
Ben Abbott <bpabbott@mac.com>
parents: 13074
diff changeset
607 %! assert (q * a * z, aa, norm (aa) * 1e-14);
4ffea87ad71b codesprint: Fix tolerance for qz.cc tests.
Ben Abbott <bpabbott@mac.com>
parents: 13074
diff changeset
608 %! assert (q * b * z, bb, norm (bb) * 1e-14);
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
609 %! assert (a * v * diag (diag (bb)), b * v * diag (diag (aa)), 1e-13);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
610 %! assert (diag (diag (bb)) * w' * a, diag (diag (aa)) * w' * b, 1e-13);
13074
0b789a03bde1 codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12642
diff changeset
611
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
612 ## An example with real matrices where some generalized eigenvalues are
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
613 ## complex, so the real and complex QZ differ.
25570
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
614 %!test
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
615 %! A = [ -1.03428 0.24929 0.43205 -0.12860;
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
616 %! 1.16228 0.27870 2.12954 0.69250;
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
617 %! -0.51524 -0.34939 -0.77820 2.13721;
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
618 %! -1.32941 2.11870 0.72005 1.00835 ];
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
619 %! B = [ 1.407302 -0.632956 -0.360628 0.068534;
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
620 %! 0.149898 0.298248 0.991777 0.023652;
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
621 %! 0.169281 -0.405205 -1.775834 1.511730;
da2b38b7a077 Use LAPACK DTGSEN subrouting for ordered-qz computation (bug #53761)
Sébastien Villemot <sebastien@debian.org>
parents: 25106
diff changeset
622 %! 0.717770 1.291390 -1.766607 -0.531352 ];
31976
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
623 %! [CAA, CBB, CQ, CZ, CV, CW] = qz (A, B, 'complex');
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
624 %! assert (iscomplex (CAA) && iscomplex (CBB) && iscomplex (CQ) ...
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
625 %! && iscomplex (CZ) && iscomplex (CV) && iscomplex (CW));
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
626 %! assert (istriu (CAA) && istriu (CBB));
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
627 %! assert (CQ * CQ', eye(4), 1e-14);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
628 %! assert (CZ * CZ', eye(4), 1e-14);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
629 %! assert (CQ * A * CZ, CAA, norm (CAA) * 1e-14);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
630 %! assert (CQ * B * CZ, CBB, norm (CBB) * 1e-14);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
631 %! assert (A * CV * diag (diag (CBB)), B * CV * diag (diag (CAA)), 1e-13);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
632 %! assert (diag (diag (CBB)) * CW' * A, diag (diag (CAA)) * CW' * B, 1e-13);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
633 %! [AA, BB, Q, Z, V, W] = qz (A, B, 'real');
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
634 %! assert (isreal (AA) && isreal (BB) && isreal (Q) && isreal (Z) ...
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
635 %! && isreal (V) && isreal (W));
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
636 %! assert (istriu (BB));
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
637 %! assert (Q * Q', eye(4), 1e-14);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
638 %! assert (Z * Z', eye(4), 1e-14);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
639 %! assert (Q * A * Z, AA, norm (AA) * 1e-14);
e6b24d0d485c qz: Handle input and output arguments in a Matlab-compatible way.
Sébastien Villemot <sebastien@debian.org>
parents: 31827
diff changeset
640 %! assert (Q * B * Z, BB, norm (BB) * 1e-14);
31977
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
641
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
642 ## Test input validation
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
643 %!error <Invalid call> qz ()
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
644 %!error <Invalid call> qz (1)
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
645 %!error <Invalid call> qz (1,2,3,4)
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
646 %!error <Invalid call> [r1,r2,r3,r4,r5,r6,r7] = qz (1,2)
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
647 %!error <OPT must be a string> qz (1,2, 3)
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
648 %!error <OPT must be 'real' or 'complex'> qz (1,2, 'foobar')
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
649 %!warning <ignoring 'real' option with complex matrices> qz (2i, 3, 'real');
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
650 %!warning <ignoring 'real' option with complex matrices> qz (2, 3i, 'real');
f5c0a0754da1 qz: code cleanup for cset e6b24d0d485c.
Rik <rik@octave.org>
parents: 31976
diff changeset
651
13074
0b789a03bde1 codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12642
diff changeset
652 */
29958
32c3a5805893 move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents: 29558
diff changeset
653
31605
e88a07dec498 maint: Use macros to begin/end C++ namespaces.
Rik <rik@octave.org>
parents: 30564
diff changeset
654 OCTAVE_END_NAMESPACE(octave)