Mercurial > octave
annotate libinterp/corefcn/qz.cc @ 20853:1142cf6abc0d
2015 Code Sprint: remove class of function from docstring for all C++ files.
author | Rik <rik@octave.org> |
---|---|
date | Sat, 12 Dec 2015 07:40:03 -0800 |
parents | 35241c4b696c |
children | c07bee629973 |
rev | line source |
---|---|
3183 | 1 /* |
2 | |
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19403
diff
changeset
|
3 Copyright (C) 1998-2015 A. S. Hodel |
3183 | 4 |
5 This file is part of Octave. | |
6 | |
7 Octave is free software; you can redistribute it and/or modify it | |
8 under the terms of the GNU General Public License as published by the | |
7016 | 9 Free Software Foundation; either version 3 of the License, or (at your |
10 option) any later version. | |
3183 | 11 |
12 Octave is distributed in the hope that it will be useful, but WITHOUT | |
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
15 for more details. | |
16 | |
17 You should have received a copy of the GNU General Public License | |
7016 | 18 along with Octave; see the file COPYING. If not, see |
19 <http://www.gnu.org/licenses/>. | |
3183 | 20 |
21 */ | |
22 | |
23 // Generalized eigenvalue balancing via LAPACK | |
3911 | 24 |
25 // Author: A. S. Hodel <scotte@eng.auburn.edu> | |
3183 | 26 |
27 #undef DEBUG | |
28 #undef DEBUG_SORT | |
29 #undef DEBUG_EIG | |
30 | |
9786
2c279308f6ab
fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
31 #ifdef HAVE_CONFIG_H |
2c279308f6ab
fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
32 #include <config.h> |
2c279308f6ab
fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
33 #endif |
3183 | 34 |
35 #include <cfloat> | |
4051 | 36 |
3523 | 37 #include <iostream> |
4051 | 38 #include <iomanip> |
3183 | 39 |
40 #include "CmplxQRP.h" | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
41 #include "CmplxQR.h" |
3183 | 42 #include "dbleQR.h" |
4153 | 43 #include "f77-fcn.h" |
7231 | 44 #include "lo-math.h" |
4153 | 45 #include "quit.h" |
46 | |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14876
diff
changeset
|
47 #include "defun.h" |
3183 | 48 #include "error.h" |
49 #include "gripes.h" | |
50 #include "oct-obj.h" | |
51 #include "oct-map.h" | |
52 #include "ov.h" | |
53 #include "pager.h" | |
3185 | 54 #if defined (DEBUG) || defined (DEBUG_SORT) |
3183 | 55 #include "pr-output.h" |
56 #endif | |
57 #include "symtab.h" | |
58 #include "utils.h" | |
59 #include "variables.h" | |
60 | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
61 typedef octave_idx_type (*sort_function) (const octave_idx_type& LSIZE, |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
62 const double& ALPHA, |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
63 const double& BETA, const double& S, |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
64 const double& P); |
3183 | 65 |
66 extern "C" | |
67 { | |
4552 | 68 F77_RET_T |
69 F77_FUNC (dggbal, DGGBAL) (F77_CONST_CHAR_ARG_DECL, | |
10311 | 70 const octave_idx_type& N, double* A, |
71 const octave_idx_type& LDA, double* B, | |
72 const octave_idx_type& LDB, octave_idx_type& ILO, | |
73 octave_idx_type& IHI, double* LSCALE, | |
74 double* RSCALE, double* WORK, | |
75 octave_idx_type& INFO | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
76 F77_CHAR_ARG_LEN_DECL); |
3183 | 77 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
78 F77_RET_T |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
79 F77_FUNC (zggbal, ZGGBAL) (F77_CONST_CHAR_ARG_DECL, |
10311 | 80 const octave_idx_type& N, Complex* A, |
81 const octave_idx_type& LDA, Complex* B, | |
82 const octave_idx_type& LDB, octave_idx_type& ILO, | |
83 octave_idx_type& IHI, double* LSCALE, | |
84 double* RSCALE, double* WORK, | |
85 octave_idx_type& INFO | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
86 F77_CHAR_ARG_LEN_DECL); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
87 |
4552 | 88 F77_RET_T |
89 F77_FUNC (dggbak, DGGBAK) (F77_CONST_CHAR_ARG_DECL, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
90 F77_CONST_CHAR_ARG_DECL, |
10311 | 91 const octave_idx_type& N, |
92 const octave_idx_type& ILO, | |
93 const octave_idx_type& IHI, | |
94 const double* LSCALE, const double* RSCALE, | |
95 octave_idx_type& M, double* V, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
96 const octave_idx_type& LDV, octave_idx_type& INFO |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
97 F77_CHAR_ARG_LEN_DECL |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
98 F77_CHAR_ARG_LEN_DECL); |
3183 | 99 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
100 F77_RET_T |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
101 F77_FUNC (zggbak, ZGGBAK) (F77_CONST_CHAR_ARG_DECL, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
102 F77_CONST_CHAR_ARG_DECL, |
10311 | 103 const octave_idx_type& N, |
104 const octave_idx_type& ILO, | |
105 const octave_idx_type& IHI, | |
106 const double* LSCALE, const double* RSCALE, | |
107 octave_idx_type& M, Complex* V, | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
108 const octave_idx_type& LDV, octave_idx_type& INFO |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
109 F77_CHAR_ARG_LEN_DECL |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
110 F77_CHAR_ARG_LEN_DECL); |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
111 |
4552 | 112 F77_RET_T |
113 F77_FUNC (dgghrd, DGGHRD) (F77_CONST_CHAR_ARG_DECL, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
114 F77_CONST_CHAR_ARG_DECL, |
10311 | 115 const octave_idx_type& N, |
116 const octave_idx_type& ILO, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
117 const octave_idx_type& IHI, double* A, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
118 const octave_idx_type& LDA, double* B, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
119 const octave_idx_type& LDB, double* Q, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
120 const octave_idx_type& LDQ, double* Z, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
121 const octave_idx_type& LDZ, octave_idx_type& INFO |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
122 F77_CHAR_ARG_LEN_DECL |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
123 F77_CHAR_ARG_LEN_DECL); |
3183 | 124 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
125 F77_RET_T |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
126 F77_FUNC (zgghrd, ZGGHRD) (F77_CONST_CHAR_ARG_DECL, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
127 F77_CONST_CHAR_ARG_DECL, |
10311 | 128 const octave_idx_type& N, |
129 const octave_idx_type& ILO, | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
130 const octave_idx_type& IHI, Complex* A, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
131 const octave_idx_type& LDA, Complex* B, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
132 const octave_idx_type& LDB, Complex* Q, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
133 const octave_idx_type& LDQ, Complex* Z, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
134 const octave_idx_type& LDZ, octave_idx_type& INFO |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
135 F77_CHAR_ARG_LEN_DECL |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
136 F77_CHAR_ARG_LEN_DECL); |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
137 |
4552 | 138 F77_RET_T |
139 F77_FUNC (dhgeqz, DHGEQZ) (F77_CONST_CHAR_ARG_DECL, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
140 F77_CONST_CHAR_ARG_DECL, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
141 F77_CONST_CHAR_ARG_DECL, |
10311 | 142 const octave_idx_type& N, |
143 const octave_idx_type& ILO, | |
144 const octave_idx_type& IHI, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
145 double* A, const octave_idx_type& LDA, double* B, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
146 const octave_idx_type& LDB, double* ALPHAR, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
147 double* ALPHAI, double* BETA, double* Q, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
148 const octave_idx_type& LDQ, double* Z, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
149 const octave_idx_type& LDZ, double* WORK, |
10311 | 150 const octave_idx_type& LWORK, |
151 octave_idx_type& INFO | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
152 F77_CHAR_ARG_LEN_DECL |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
153 F77_CHAR_ARG_LEN_DECL |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
154 F77_CHAR_ARG_LEN_DECL); |
3183 | 155 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
156 F77_RET_T |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
157 F77_FUNC (zhgeqz, ZHGEQZ) (F77_CONST_CHAR_ARG_DECL, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
158 F77_CONST_CHAR_ARG_DECL, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
159 F77_CONST_CHAR_ARG_DECL, |
10311 | 160 const octave_idx_type& N, |
161 const octave_idx_type& ILO, | |
162 const octave_idx_type& IHI, | |
163 Complex* A, const octave_idx_type& LDA, | |
164 Complex* B, const octave_idx_type& LDB, | |
165 Complex* ALPHA, Complex* BETA, Complex* CQ, | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
166 const octave_idx_type& LDQ, |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
167 Complex* CZ, const octave_idx_type& LDZ, |
10311 | 168 Complex* WORK, const octave_idx_type& LWORK, |
169 double* RWORK, octave_idx_type& INFO | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
170 F77_CHAR_ARG_LEN_DECL |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
171 F77_CHAR_ARG_LEN_DECL |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
172 F77_CHAR_ARG_LEN_DECL); |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
173 |
4552 | 174 F77_RET_T |
10311 | 175 F77_FUNC (dlag2, DLAG2) (const double* A, const octave_idx_type& LDA, |
176 const double* B, const octave_idx_type& LDB, | |
177 const double& SAFMIN, double& SCALE1, | |
178 double& SCALE2, double& WR1, double& WR2, | |
179 double& WI); | |
3183 | 180 |
181 // Van Dooren's code (netlib.org: toms/590) for reordering | |
182 // GEP. Only processes Z, not Q. | |
4552 | 183 F77_RET_T |
10311 | 184 F77_FUNC (dsubsp, DSUBSP) (const octave_idx_type& NMAX, |
185 const octave_idx_type& N, double* A, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
186 double* B, double* Z, sort_function, |
10311 | 187 const double& EPS, octave_idx_type& NDIM, |
188 octave_idx_type& FAIL, octave_idx_type* IND); | |
3183 | 189 |
10311 | 190 // Documentation for DTGEVC incorrectly states that VR, VL are |
3183 | 191 // complex*16; they are declared in DTGEVC as double precision |
10311 | 192 // (probably a cut and paste problem fro ZTGEVC). |
4552 | 193 F77_RET_T |
194 F77_FUNC (dtgevc, DTGEVC) (F77_CONST_CHAR_ARG_DECL, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
195 F77_CONST_CHAR_ARG_DECL, |
10311 | 196 octave_idx_type* SELECT, |
197 const octave_idx_type& N, double* A, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
198 const octave_idx_type& LDA, double* B, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
199 const octave_idx_type& LDB, double* VL, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
200 const octave_idx_type& LDVL, double* VR, |
10311 | 201 const octave_idx_type& LDVR, |
202 const octave_idx_type& MM, octave_idx_type& M, | |
203 double* WORK, octave_idx_type& INFO | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
204 F77_CHAR_ARG_LEN_DECL |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
205 F77_CHAR_ARG_LEN_DECL); |
3183 | 206 |
19861
19755f4fc851
maint: Cleanup C++ code to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
207 F77_RET_T |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
208 F77_FUNC (ztgevc, ZTGEVC) (F77_CONST_CHAR_ARG_DECL, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
209 F77_CONST_CHAR_ARG_DECL, |
10311 | 210 octave_idx_type* SELECT, |
211 const octave_idx_type& N, const Complex* A, | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
212 const octave_idx_type& LDA,const Complex* B, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
213 const octave_idx_type& LDB, Complex* xVL, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
214 const octave_idx_type& LDVL, Complex* xVR, |
10311 | 215 const octave_idx_type& LDVR, |
216 const octave_idx_type& MM, octave_idx_type& M, | |
217 Complex* CWORK, double* RWORK, | |
218 octave_idx_type& INFO | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
219 F77_CHAR_ARG_LEN_DECL |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
220 F77_CHAR_ARG_LEN_DECL); |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
221 |
4552 | 222 F77_RET_T |
223 F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG_DECL, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
224 double& retval |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
225 F77_CHAR_ARG_LEN_DECL); |
3185 | 226 |
4552 | 227 F77_RET_T |
228 F77_FUNC (xdlange, XDLANGE) (F77_CONST_CHAR_ARG_DECL, | |
10311 | 229 const octave_idx_type&, |
230 const octave_idx_type&, const double*, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
231 const octave_idx_type&, double*, double& |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
232 F77_CHAR_ARG_LEN_DECL); |
3183 | 233 } |
234 | |
235 // fcrhp, fin, fout, folhp: | |
236 // routines for ordering of generalized eigenvalues | |
237 // return 1 if test is passed, 0 otherwise | |
238 // fin: |lambda| < 1 | |
239 // fout: |lambda| >= 1 | |
240 // fcrhp: real(lambda) >= 0 | |
241 // folhp: real(lambda) < 0 | |
242 | |
5275 | 243 static octave_idx_type |
244 fcrhp (const octave_idx_type& lsize, const double& alpha, | |
3185 | 245 const double& beta, const double& s, const double&) |
3183 | 246 { |
3185 | 247 if (lsize == 1) |
10311 | 248 return (alpha * beta >= 0 ? 1 : -1); |
3185 | 249 else |
3183 | 250 return (s >= 0 ? 1 : -1); |
251 } | |
3185 | 252 |
5275 | 253 static octave_idx_type |
254 fin (const octave_idx_type& lsize, const double& alpha, | |
3185 | 255 const double& beta, const double&, const double& p) |
3183 | 256 { |
5275 | 257 octave_idx_type retval; |
3183 | 258 |
3185 | 259 if (lsize == 1) |
260 retval = (fabs (alpha) < fabs (beta) ? 1 : -1); | |
261 else | |
262 retval = (fabs (p) < 1 ? 1 : -1); | |
263 | |
264 #ifdef DEBUG | |
3538 | 265 std::cout << "qz: fin: retval=" << retval << std::endl; |
3185 | 266 #endif |
267 | |
3183 | 268 return retval; |
269 } | |
3185 | 270 |
5275 | 271 static octave_idx_type |
272 folhp (const octave_idx_type& lsize, const double& alpha, | |
3185 | 273 const double& beta, const double& s, const double&) |
3183 | 274 { |
3185 | 275 if (lsize == 1) |
10311 | 276 return (alpha * beta < 0 ? 1 : -1); |
3185 | 277 else |
3183 | 278 return (s < 0 ? 1 : -1); |
279 } | |
3185 | 280 |
5275 | 281 static octave_idx_type |
282 fout (const octave_idx_type& lsize, const double& alpha, | |
3185 | 283 const double& beta, const double&, const double& p) |
3183 | 284 { |
3185 | 285 if (lsize == 1) |
286 return (fabs (alpha) >= fabs (beta) ? 1 : -1); | |
287 else | |
288 return (fabs (p) >= 1 ? 1 : -1); | |
3183 | 289 } |
290 | |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
291 |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
292 //FIXME: Matlab does not produce lambda as the first output argument. |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
293 // Compatibility problem? |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14876
diff
changeset
|
294 DEFUN (qz, args, nargout, |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
295 "-*- texinfo -*-\n\ |
20853
1142cf6abc0d
2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents:
20831
diff
changeset
|
296 @deftypefn {} {@var{lambda} =} qz (@var{A}, @var{B})\n\ |
1142cf6abc0d
2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents:
20831
diff
changeset
|
297 @deftypefnx {} {@var{lambda} =} qz (@var{A}, @var{B}, @var{opt})\n\ |
11593
1577c6f80926
Use non-breaking spaces between certain adjectives and their nouns in docstrings.
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
298 QZ@tie{}decomposition of the generalized eigenvalue problem\n\ |
20172
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
299 (@math{A x = s B x}).\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
300 \n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
301 There are three ways to call this function:\n\ |
3372 | 302 @enumerate\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
303 @item @code{@var{lambda} = qz (@var{A}, @var{B})}\n\ |
3185 | 304 \n\ |
5016 | 305 Computes the generalized eigenvalues\n\ |
306 @tex\n\ | |
307 $\\lambda$\n\ | |
308 @end tex\n\ | |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
309 @ifnottex\n\ |
5016 | 310 @var{lambda}\n\ |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
311 @end ifnottex\n\ |
5016 | 312 of @math{(A - s B)}.\n\ |
10840 | 313 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
314 @item @code{[AA, BB, Q, Z, V, W, @var{lambda}] = qz (@var{A}, @var{B})}\n\ |
3185 | 315 \n\ |
20172
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
316 Computes QZ@tie{}decomposition, generalized eigenvectors, and generalized\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
317 eigenvalues of @math{(A - s B)}\n\ |
5016 | 318 @tex\n\ |
319 $$ AV = BV{ \\rm diag }(\\lambda) $$\n\ | |
320 $$ W^T A = { \\rm diag }(\\lambda)W^T B $$\n\ | |
321 $$ AA = Q^T AZ, BB = Q^T BZ $$\n\ | |
322 @end tex\n\ | |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
323 @ifnottex\n\ |
10840 | 324 \n\ |
3372 | 325 @example\n\ |
326 @group\n\ | |
5481 | 327 \n\ |
14360
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
328 A * V = B * V * diag (@var{lambda})\n\ |
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
329 W' * A = diag (@var{lambda}) * W' * B\n\ |
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
330 AA = Q * A * Z, BB = Q * B * Z\n\ |
5481 | 331 \n\ |
3372 | 332 @end group\n\ |
333 @end example\n\ | |
10840 | 334 \n\ |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
335 @end ifnottex\n\ |
5016 | 336 with @var{Q} and @var{Z} orthogonal (unitary)= @var{I}\n\ |
3185 | 337 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
338 @item @code{[AA,BB,Z@{, @var{lambda}@}] = qz (@var{A}, @var{B}, @var{opt})}\n\ |
3185 | 339 \n\ |
20172
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
340 As in form [2], but allows ordering of generalized eigenpairs for, e.g.,\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
341 solution of discrete time algebraic Riccati equations. Form 3 is not\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
342 available for complex matrices, and does not compute the generalized\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
343 eigenvectors @var{V}, @var{W}, nor the orthogonal matrix @var{Q}.\n\ |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10350
diff
changeset
|
344 \n\ |
3372 | 345 @table @var\n\ |
346 @item opt\n\ | |
20172
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
347 for ordering eigenvalues of the @nospell{GEP} pencil. The leading block of\n\ |
4f45eaf83908
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19861
diff
changeset
|
348 the revised pencil contains all eigenvalues that satisfy:\n\ |
14360
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
349 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
350 @table @asis\n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
351 @item @qcode{\"N\"}\n\ |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
11593
diff
changeset
|
352 = unordered (default)\n\ |
3183 | 353 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
354 @item @qcode{\"S\"}\n\ |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
11593
diff
changeset
|
355 = small: leading block has all |lambda| @leq{} 1\n\ |
3185 | 356 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
357 @item @qcode{\"B\"}\n\ |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
11593
diff
changeset
|
358 = big: leading block has all |lambda| @geq{} 1\n\ |
3372 | 359 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
360 @item @qcode{\"-\"}\n\ |
5481 | 361 = negative real part: leading block has all eigenvalues\n\ |
362 in the open left half-plane\n\ | |
3372 | 363 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
364 @item @qcode{\"+\"}\n\ |
9209
923c7cb7f13f
Simplify TeXinfo files by eliminating redundant @iftex followed by @tex construction.
Rik <rdrider0-list@yahoo.com>
parents:
9064
diff
changeset
|
365 = non-negative real part: leading block has all eigenvalues\n\ |
5481 | 366 in the closed right half-plane\n\ |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8927
diff
changeset
|
367 @end table\n\ |
3372 | 368 @end table\n\ |
369 @end enumerate\n\ | |
3183 | 370 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
371 Note: @code{qz} performs permutation balancing, but not scaling\n\ |
17097
e7a059a9a644
doc: Use XREF as anchor prefix in documentation for clearer results in Info viewer.
Rik <rik@octave.org>
parents:
16920
diff
changeset
|
372 (@pxref{XREFbalance}). The order of output arguments was selected for\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
373 compatibility with @sc{matlab}.\n\ |
16920
53eaa83e4181
doc: Add seealso links between various factorization forms.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
374 @seealso{eig, balance, lu, chol, hess, qr, qzhess, schur, svd}\n\ |
3372 | 375 @end deftypefn") |
3183 | 376 { |
377 octave_value_list retval; | |
378 int nargin = args.length (); | |
379 | |
3185 | 380 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
381 std::cout << "qz: nargin = " << nargin |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
382 << ", nargout = " << nargout << std::endl; |
3185 | 383 #endif |
3183 | 384 |
3185 | 385 if (nargin < 2 || nargin > 3 || nargout > 7) |
20802
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20700
diff
changeset
|
386 print_usage (); |
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20700
diff
changeset
|
387 |
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20700
diff
changeset
|
388 if (nargin == 3 && (nargout < 3 || nargout > 4)) |
8bb38ba1bad6
eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents:
20700
diff
changeset
|
389 error ("qz: invalid number of output arguments for form [3] call"); |
3183 | 390 |
3185 | 391 #ifdef DEBUG |
3538 | 392 std::cout << "qz: determine ordering option" << std::endl; |
3185 | 393 #endif |
3183 | 394 |
10311 | 395 // Determine ordering option. |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
396 volatile char ord_job = 0; |
3183 | 397 static double safmin; |
3185 | 398 |
399 if (nargin == 2) | |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
400 ord_job = 'N'; |
3185 | 401 else |
402 { | |
20700
68e3a747ca02
rename octave_value value extractors that accept error message args
John W. Eaton <jwe@octave.org>
parents:
20581
diff
changeset
|
403 std::string tmp = 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
|
404 |
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
405 if (! tmp.empty ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
406 ord_job = tmp[0]; |
3183 | 407 |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
408 if (! (ord_job == 'N' || ord_job == 'n' |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
409 || ord_job == 'S' || ord_job == 's' |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
410 || ord_job == 'B' || ord_job == 'b' |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
411 || ord_job == '+' || ord_job == '-')) |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
412 error ("qz: invalid order option"); |
3185 | 413 |
414 // overflow constant required by dlag2 | |
4552 | 415 F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG2 ("S", 1), |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
416 safmin |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
417 F77_CHAR_ARG_LEN (1)); |
3183 | 418 |
3185 | 419 #ifdef DEBUG_EIG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
420 std::cout << "qz: initial value of safmin=" |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
421 << setiosflags (std::ios::scientific) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
422 << safmin << std::endl; |
3185 | 423 #endif |
3183 | 424 |
10311 | 425 // Some machines (e.g., DEC alpha) get safmin = 0; |
426 // for these, use eps instead to avoid problems in dlag2. | |
3185 | 427 if (safmin == 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
428 { |
3185 | 429 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
430 std::cout << "qz: DANGER WILL ROBINSON: safmin is 0!" << std::endl; |
3185 | 431 #endif |
3183 | 432 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
433 F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG2 ("E", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
434 safmin |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
435 F77_CHAR_ARG_LEN (1)); |
3185 | 436 |
437 #ifdef DEBUG_EIG | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
438 std::cout << "qz: safmin set to " |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
439 << setiosflags (std::ios::scientific) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
440 << safmin << std::endl; |
3185 | 441 #endif |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
442 } |
3183 | 443 } |
444 | |
3185 | 445 #ifdef DEBUG |
3538 | 446 std::cout << "qz: check argument 1" << std::endl; |
3185 | 447 #endif |
3183 | 448 |
10311 | 449 // Argument 1: check if it's o.k. dimensioned. |
5275 | 450 octave_idx_type nn = args(0).rows (); |
3185 | 451 |
452 #ifdef DEBUG | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
453 std::cout << "argument 1 dimensions: (" |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
454 << nn << "," << args(0).columns () << ")" |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
455 << std::endl; |
3185 | 456 #endif |
457 | |
458 int arg_is_empty = empty_arg ("qz", nn, args(0).columns ()); | |
459 | |
3183 | 460 if (arg_is_empty < 0) |
3185 | 461 { |
462 gripe_empty_arg ("qz: parameter 1", 0); | |
463 return retval; | |
464 } | |
3183 | 465 else if (arg_is_empty > 0) |
3185 | 466 { |
467 gripe_empty_arg ("qz: parameter 1; continuing", 0); | |
468 return octave_value_list (2, Matrix ()); | |
469 } | |
470 else if (args(0).columns () != nn) | |
471 { | |
472 gripe_square_matrix_required ("qz"); | |
473 return retval; | |
474 } | |
3183 | 475 |
10311 | 476 // Argument 1: dimensions look good; get the value. |
3183 | 477 Matrix aa; |
478 ComplexMatrix caa; | |
3185 | 479 |
480 if (args(0).is_complex_type ()) | |
3183 | 481 caa = args(0).complex_matrix_value (); |
3185 | 482 else |
3183 | 483 aa = args(0).matrix_value (); |
3185 | 484 |
485 #ifdef DEBUG | |
3538 | 486 std::cout << "qz: check argument 2" << std::endl; |
3185 | 487 #endif |
3183 | 488 |
10311 | 489 // Extract argument 2 (bb, or cbb if complex). |
3185 | 490 if ((nn != args(1).columns ()) || (nn != args(1).rows ())) |
491 { | |
492 gripe_nonconformant (); | |
493 return retval; | |
494 } | |
495 | |
3183 | 496 Matrix bb; |
497 ComplexMatrix cbb; | |
3185 | 498 |
499 if (args(1).is_complex_type ()) | |
3183 | 500 cbb = args(1).complex_matrix_value (); |
501 else | |
502 bb = args(1).matrix_value (); | |
3185 | 503 |
3183 | 504 // Both matrices loaded, now let's check what kind of arithmetic: |
10311 | 505 // declared volatile to avoid compiler warnings about long jumps, |
506 // vforks. | |
3185 | 507 |
10311 | 508 volatile int complex_case |
3185 | 509 = (args(0).is_complex_type () || args(1).is_complex_type ()); |
10311 | 510 |
3185 | 511 if (nargin == 3 && complex_case) |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
512 error ("qz: cannot re-order complex qz decomposition"); |
3183 | 513 |
10311 | 514 // First, declare variables used in both the real and complex case. |
3183 | 515 Matrix QQ(nn,nn), ZZ(nn,nn), VR(nn,nn), VL(nn,nn); |
516 RowVector alphar(nn), alphai(nn), betar(nn); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
517 ComplexRowVector xalpha(nn), xbeta(nn); |
3185 | 518 ComplexMatrix CQ(nn,nn), CZ(nn,nn), CVR(nn,nn), CVL(nn,nn); |
5275 | 519 octave_idx_type ilo, ihi, info; |
3185 | 520 char compq = (nargout >= 3 ? 'V' : 'N'); |
14876
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
521 char compz = ((nargout >= 4 || nargin == 3)? 'V' : 'N'); |
3183 | 522 |
10311 | 523 // Initialize Q, Z to identity if we need either of them. |
3185 | 524 if (compq == 'V' || compz == 'V') |
5275 | 525 for (octave_idx_type ii = 0; ii < nn; ii++) |
526 for (octave_idx_type jj = 0; jj < nn; jj++) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
527 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
528 OCTAVE_QUIT; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
529 QQ(ii,jj) = ZZ(ii,jj) = (ii == jj ? 1.0 : 0.0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
530 } |
3183 | 531 |
10311 | 532 // Always perform permutation balancing. |
4552 | 533 const char bal_job = 'P'; |
10311 | 534 RowVector lscale (nn), rscale (nn), work (6 * nn), rwork (nn); |
3183 | 535 |
3185 | 536 if (complex_case) |
537 { | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
538 #ifdef DEBUG |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
539 if (compq == 'V') |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
540 std::cout << "qz: performing balancing; CQ=" << std::endl |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
541 << CQ << std::endl; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
542 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
543 if (args(0).is_real_type ()) |
10311 | 544 caa = ComplexMatrix (aa); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
545 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
546 if (args(1).is_real_type ()) |
10311 | 547 cbb = ComplexMatrix (bb); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
548 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
549 if (compq == 'V') |
10311 | 550 CQ = ComplexMatrix (QQ); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
551 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
552 if (compz == 'V') |
10311 | 553 CZ = ComplexMatrix (ZZ); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
554 |
10311 | 555 F77_XFCN (zggbal, ZGGBAL, |
556 (F77_CONST_CHAR_ARG2 (&bal_job, 1), | |
557 nn, caa.fortran_vec (), nn, cbb.fortran_vec (), | |
558 nn, ilo, ihi, lscale.fortran_vec (), | |
559 rscale.fortran_vec (), work.fortran_vec (), info | |
560 F77_CHAR_ARG_LEN (1))); | |
3185 | 561 } |
3183 | 562 else |
3185 | 563 { |
564 #ifdef DEBUG | |
565 if (compq == 'V') | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
566 std::cout << "qz: performing balancing; QQ=" << std::endl |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
567 << QQ << std::endl; |
3185 | 568 #endif |
3183 | 569 |
3185 | 570 F77_XFCN (dggbal, DGGBAL, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
571 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
572 nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
573 nn, ilo, ihi, lscale.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
574 rscale.fortran_vec (), work.fortran_vec (), info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
575 F77_CHAR_ARG_LEN (1))); |
3185 | 576 } |
3183 | 577 |
578 // Since we just want the balancing matrices, we can use dggbal | |
10311 | 579 // for both the real and complex cases; left first |
3185 | 580 |
10311 | 581 #if 0 |
582 if (compq == 'V') | |
3185 | 583 { |
584 F77_XFCN (dggbak, DGGBAK, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
585 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
586 F77_CONST_CHAR_ARG2 ("L", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
587 nn, ilo, ihi, lscale.data (), rscale.data (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
588 nn, QQ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
589 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
590 F77_CHAR_ARG_LEN (1))); |
3183 | 591 |
3185 | 592 #ifdef DEBUG |
593 if (compq == 'V') | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
594 std::cout << "qz: balancing done; QQ=" << std::endl << QQ << std::endl; |
3185 | 595 #endif |
3183 | 596 } |
597 | |
598 // then right | |
3185 | 599 if (compz == 'V') |
600 { | |
4552 | 601 F77_XFCN (dggbak, DGGBAK, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
602 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
603 F77_CONST_CHAR_ARG2 ("R", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
604 nn, ilo, ihi, lscale.data (), rscale.data (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
605 nn, ZZ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
606 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
607 F77_CHAR_ARG_LEN (1))); |
3183 | 608 |
3185 | 609 #ifdef DEBUG |
610 if (compz == 'V') | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
611 std::cout << "qz: balancing done; ZZ=" << std::endl << ZZ << std::endl; |
3185 | 612 #endif |
10311 | 613 } |
614 #endif | |
3183 | 615 |
616 static char qz_job; | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
617 qz_job = (nargout < 2 ? 'E' : 'S'); |
3185 | 618 |
3183 | 619 if (complex_case) |
3185 | 620 { |
10311 | 621 // Complex case. |
622 | |
623 // The QR decomposition of cbb. | |
624 ComplexQR cbqr (cbb); | |
625 // The R matrix of QR decomposition for cbb. | |
626 cbb = cbqr.R (); | |
627 // (Q*)caa for following work. | |
628 caa = (cbqr.Q ().hermitian ()) * caa; | |
629 CQ = CQ * cbqr.Q (); | |
630 | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
631 F77_XFCN (zgghrd, ZGGHRD, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
632 (F77_CONST_CHAR_ARG2 (&compq, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
633 F77_CONST_CHAR_ARG2 (&compz, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
634 nn, ilo, ihi, caa.fortran_vec (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
635 nn, cbb.fortran_vec (), nn, CQ.fortran_vec (), nn, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
636 CZ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
637 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
638 F77_CHAR_ARG_LEN (1))); |
10311 | 639 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
640 ComplexRowVector cwork (1 * nn); |
10311 | 641 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
642 F77_XFCN (zhgeqz, ZHGEQZ, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
643 (F77_CONST_CHAR_ARG2 (&qz_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
644 F77_CONST_CHAR_ARG2 (&compq, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
645 F77_CONST_CHAR_ARG2 (&compz, 1), |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
646 nn, ilo, ihi, |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
647 caa.fortran_vec (), nn, |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
648 cbb.fortran_vec (),nn, |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
649 xalpha.fortran_vec (), xbeta.fortran_vec (), |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
650 CQ.fortran_vec (), nn, |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
651 CZ.fortran_vec (), nn, |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
652 cwork.fortran_vec (), nn, rwork.fortran_vec (), info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
653 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
654 F77_CHAR_ARG_LEN (1) |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
655 F77_CHAR_ARG_LEN (1))); |
3185 | 656 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
657 if (compq == 'V') |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
658 { |
10311 | 659 // Left eigenvector. |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
660 F77_XFCN (zggbak, ZGGBAK, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
661 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
662 F77_CONST_CHAR_ARG2 ("L", 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
663 nn, ilo, ihi, lscale.data (), rscale.data (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
664 nn, CQ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
665 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
666 F77_CHAR_ARG_LEN (1))); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
667 } |
3185 | 668 |
10311 | 669 // Right eigenvector. |
3185 | 670 if (compz == 'V') |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
671 { |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
672 F77_XFCN (zggbak, ZGGBAK, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
673 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
674 F77_CONST_CHAR_ARG2 ("R", 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
675 nn, ilo, ihi, lscale.data (), rscale.data (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
676 nn, CZ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
677 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
678 F77_CHAR_ARG_LEN (1))); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
679 } |
3185 | 680 |
681 } | |
10311 | 682 else |
3185 | 683 { |
684 #ifdef DEBUG | |
3538 | 685 std::cout << "qz: peforming qr decomposition of bb" << std::endl; |
3185 | 686 #endif |
3183 | 687 |
10311 | 688 // Compute the QR factorization of bb. |
3185 | 689 QR bqr (bb); |
690 | |
691 #ifdef DEBUG | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
692 std::cout << "qz: qr (bb) done; now peforming qz decomposition" |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
693 << std::endl; |
3185 | 694 #endif |
3183 | 695 |
3185 | 696 bb = bqr.R (); |
697 | |
698 #ifdef DEBUG | |
3538 | 699 std::cout << "qz: extracted bb" << std::endl; |
3185 | 700 #endif |
3183 | 701 |
10311 | 702 aa = (bqr.Q ()).transpose () * aa; |
3185 | 703 |
704 #ifdef DEBUG | |
3538 | 705 std::cout << "qz: updated aa " << std::endl; |
706 std::cout << "bqr.Q () = " << std::endl << bqr.Q () << std::endl; | |
3183 | 707 |
3185 | 708 if (compq == 'V') |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
709 std::cout << "QQ =" << QQ << std::endl; |
3185 | 710 #endif |
3183 | 711 |
3185 | 712 if (compq == 'V') |
10311 | 713 QQ = QQ * bqr.Q (); |
3183 | 714 |
3185 | 715 #ifdef DEBUG |
3538 | 716 std::cout << "qz: precursors done..." << std::endl; |
3185 | 717 #endif |
3183 | 718 |
3185 | 719 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
720 std::cout << "qz: compq = " << compq << ", compz = " << compz |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
721 << std::endl; |
3185 | 722 #endif |
3183 | 723 |
10311 | 724 // Reduce to generalized hessenberg form. |
3185 | 725 F77_XFCN (dgghrd, DGGHRD, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
726 (F77_CONST_CHAR_ARG2 (&compq, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
727 F77_CONST_CHAR_ARG2 (&compz, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
728 nn, ilo, ihi, aa.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
729 nn, bb.fortran_vec (), nn, QQ.fortran_vec (), nn, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
730 ZZ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
731 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
732 F77_CHAR_ARG_LEN (1))); |
3183 | 733 |
10311 | 734 // Check if just computing generalized eigenvalues or if we're |
735 // actually computing the decomposition. | |
3183 | 736 |
10311 | 737 // Reduce to generalized Schur form. |
3185 | 738 F77_XFCN (dhgeqz, DHGEQZ, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
739 (F77_CONST_CHAR_ARG2 (&qz_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
740 F77_CONST_CHAR_ARG2 (&compq, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
741 F77_CONST_CHAR_ARG2 (&compz, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
742 nn, ilo, ihi, aa.fortran_vec (), nn, bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
743 nn, alphar.fortran_vec (), alphai.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
744 betar.fortran_vec (), QQ.fortran_vec (), nn, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
745 ZZ.fortran_vec (), nn, work.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
746 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
747 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
748 F77_CHAR_ARG_LEN (1))); |
10311 | 749 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
750 if (compq == 'V') |
10311 | 751 { |
752 F77_XFCN (dggbak, DGGBAK, | |
753 (F77_CONST_CHAR_ARG2 (&bal_job, 1), | |
754 F77_CONST_CHAR_ARG2 ("L", 1), | |
755 nn, ilo, ihi, lscale.data (), rscale.data (), | |
756 nn, QQ.fortran_vec (), nn, info | |
757 F77_CHAR_ARG_LEN (1) | |
758 F77_CHAR_ARG_LEN (1))); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
759 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
760 #ifdef DEBUG |
10311 | 761 if (compq == 'V') |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
762 std::cout << "qz: balancing done; QQ=" << std::endl |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
763 << QQ << std::endl; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
764 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
765 } |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
766 |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
767 // then right |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
768 if (compz == 'V') |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
769 { |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
770 F77_XFCN (dggbak, DGGBAK, |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
771 (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
|
772 F77_CONST_CHAR_ARG2 ("R", 1), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
773 nn, ilo, ihi, lscale.data (), rscale.data (), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
774 nn, ZZ.fortran_vec (), nn, info |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
775 F77_CHAR_ARG_LEN (1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
776 F77_CHAR_ARG_LEN (1))); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
777 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
778 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
779 if (compz == 'V') |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
780 std::cout << "qz: balancing done; ZZ=" << std::endl |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
781 << ZZ << std::endl; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
782 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
783 } |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
784 |
3185 | 785 } |
3183 | 786 |
10311 | 787 // Order the QZ decomposition? |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
788 if (! (ord_job == 'N' || ord_job == 'n')) |
3183 | 789 { |
3185 | 790 if (complex_case) |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
791 // Probably not needed, but better be safe. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
792 error ("qz: cannot re-order complex qz decomposition"); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
793 |
3185 | 794 #ifdef DEBUG_SORT |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
795 std::cout << "qz: ordering eigenvalues: ord_job = " |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
796 << ord_job << std::endl; |
3185 | 797 #endif |
3183 | 798 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
799 // Declared static to avoid vfork/long jump compiler complaints. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
800 static sort_function sort_test; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
801 sort_test = 0; |
3183 | 802 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
803 switch (ord_job) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
804 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
805 case 'S': |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
806 case 's': |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
807 sort_test = &fin; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
808 break; |
3183 | 809 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
810 case 'B': |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
811 case 'b': |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
812 sort_test = &fout; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
813 break; |
3183 | 814 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
815 case '+': |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
816 sort_test = &fcrhp; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
817 break; |
3183 | 818 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
819 case '-': |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
820 sort_test = &folhp; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
821 break; |
3185 | 822 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
823 default: |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
824 // Invalid order option (should never happen, since we |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
825 // checked the options at the top). |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
826 panic_impossible (); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
827 break; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
828 } |
3183 | 829 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
830 octave_idx_type ndim, fail; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
831 double inf_norm; |
3185 | 832 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
833 F77_XFCN (xdlange, XDLANGE, |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
834 (F77_CONST_CHAR_ARG2 ("I", 1), |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
835 nn, nn, aa.data (), nn, work.fortran_vec (), inf_norm |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
836 F77_CHAR_ARG_LEN (1))); |
3185 | 837 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
838 double eps = std::numeric_limits<double>::epsilon () * inf_norm * nn; |
3185 | 839 |
840 #ifdef DEBUG_SORT | |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
841 std::cout << "qz: calling dsubsp: aa=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
842 octave_print_internal (std::cout, aa, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
843 std::cout << std::endl << "bb=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
844 octave_print_internal (std::cout, bb, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
845 if (compz == 'V') |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
846 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
847 std::cout << std::endl << "ZZ=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
848 octave_print_internal (std::cout, ZZ, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
849 } |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
850 std::cout << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
851 std::cout << "alphar = " << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
852 octave_print_internal (std::cout, (Matrix) alphar, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
853 std::cout << std::endl << "alphai = " << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
854 octave_print_internal (std::cout, (Matrix) alphai, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
855 std::cout << std::endl << "beta = " << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
856 octave_print_internal (std::cout, (Matrix) betar, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
857 std::cout << std::endl; |
3185 | 858 #endif |
859 | |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
860 Array<octave_idx_type> ind (dim_vector (nn, 1)); |
3550 | 861 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
862 F77_XFCN (dsubsp, DSUBSP, |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
863 (nn, nn, aa.fortran_vec (), bb.fortran_vec (), |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
864 ZZ.fortran_vec (), sort_test, eps, ndim, fail, |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
865 ind.fortran_vec ())); |
3185 | 866 |
867 #ifdef DEBUG | |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
868 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
869 octave_print_internal (std::cout, aa, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
870 std::cout << std::endl << "bb=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
871 octave_print_internal (std::cout, bb, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
872 if (compz == 'V') |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
873 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
874 std::cout << std::endl << "ZZ=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
875 octave_print_internal (std::cout, ZZ, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
876 } |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
877 std::cout << std::endl; |
3185 | 878 #endif |
879 | |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
880 // Manually update alphar, alphai, betar. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
881 static int jj; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
882 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
883 jj = 0; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
884 while (jj < nn) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
885 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
886 #ifdef DEBUG_EIG |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
887 std::cout << "computing gen eig #" << jj << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
888 #endif |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
889 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
890 // Number of zeros in this block. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
891 static int zcnt; |
3185 | 892 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
893 if (jj == (nn-1)) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
894 zcnt = 1; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
895 else if (aa(jj+1,jj) == 0) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
896 zcnt = 1; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
897 else zcnt = 2; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
898 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
899 if (zcnt == 1) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
900 { |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
901 // Real zero. |
3185 | 902 #ifdef DEBUG_EIG |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
903 std::cout << " single gen eig:" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
904 std::cout << " alphar(" << jj << ") = " << aa(jj,jj) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
905 << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
906 std::cout << " betar(" << jj << ") = " << bb(jj,jj) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
907 << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
908 std::cout << " alphai(" << jj << ") = 0" << std::endl; |
3185 | 909 #endif |
910 | |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
911 alphar(jj) = aa(jj,jj); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
912 alphai(jj) = 0; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
913 betar(jj) = bb(jj,jj); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
914 } |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
915 else |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
916 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
917 // Complex conjugate pair. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
918 #ifdef DEBUG_EIG |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
919 std::cout << "qz: calling dlag2:" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
920 std::cout << "safmin=" |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
921 << setiosflags (std::ios::scientific) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
922 << safmin << std::endl; |
3185 | 923 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
924 for (int idr = jj; idr <= jj+1; idr++) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
925 { |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
926 for (int idc = jj; idc <= jj+1; idc++) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
927 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
928 std::cout << "aa(" << idr << "," << idc << ")=" |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
929 << aa(idr,idc) << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
930 std::cout << "bb(" << idr << "," << idc << ")=" |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
931 << bb(idr,idc) << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
932 } |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
933 } |
3185 | 934 #endif |
935 | |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
936 // FIXME: probably should be using |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
937 // fortran_vec instead of &aa(jj,jj) here. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
938 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
939 double scale1, scale2, wr1, wr2, wi; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
940 const double *aa_ptr = aa.data () + jj * nn + jj; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
941 const double *bb_ptr = bb.data () + jj * nn + jj; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
942 F77_XFCN (dlag2, DLAG2, |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
943 (aa_ptr, nn, bb_ptr, nn, safmin, |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
944 scale1, scale2, wr1, wr2, wi)); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
945 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
946 #ifdef DEBUG_EIG |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
947 std::cout << "dlag2 returns: scale1=" << scale1 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
948 << "\tscale2=" << scale2 << std::endl |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
949 << "\twr1=" << wr1 << "\twr2=" << wr2 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
950 << "\twi=" << wi << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
951 #endif |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
952 |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
953 // Just to be safe, check if it's a real pair. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
954 if (wi == 0) |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
955 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
956 alphar(jj) = wr1; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
957 alphai(jj) = 0; |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
958 betar(jj) = scale1; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
959 alphar(jj+1) = wr2; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
960 alphai(jj+1) = 0; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
961 betar(jj+1) = scale2; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
962 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
963 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
964 { |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
965 alphar(jj) = alphar(jj+1) = wr1; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
966 alphai(jj) = -(alphai(jj+1) = wi); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
967 betar(jj) = betar(jj+1) = scale1; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
968 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
969 } |
3185 | 970 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
971 // Advance past this block. |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
972 jj += zcnt; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
973 } |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
974 |
3185 | 975 #ifdef DEBUG_SORT |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
976 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
977 octave_print_internal (std::cout, aa, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
978 std::cout << std::endl << "bb=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
979 octave_print_internal (std::cout, bb, 0); |
3185 | 980 |
20831
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
981 if (compz == 'V') |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
982 { |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
983 std::cout << std::endl << "ZZ=" << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
984 octave_print_internal (std::cout, ZZ, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
985 } |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
986 std::cout << std::endl << "qz: ndim=" << ndim << std::endl |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
987 << "fail=" << fail << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
988 std::cout << "alphar = " << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
989 octave_print_internal (std::cout, (Matrix) alphar, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
990 std::cout << std::endl << "alphai = " << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
991 octave_print_internal (std::cout, (Matrix) alphai, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
992 std::cout << std::endl << "beta = " << std::endl; |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
993 octave_print_internal (std::cout, (Matrix) betar, 0); |
35241c4b696c
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20802
diff
changeset
|
994 std::cout << std::endl; |
3185 | 995 #endif |
3183 | 996 } |
3185 | 997 |
10311 | 998 // Compute generalized eigenvalues? |
3183 | 999 ComplexColumnVector gev; |
3185 | 1000 |
1001 if (nargout < 2 || nargout == 7 || (nargin == 3 && nargout == 4)) | |
3183 | 1002 { |
3185 | 1003 if (complex_case) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1004 { |
10311 | 1005 int cnt = 0; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1006 |
10311 | 1007 for (int ii = 0; ii < nn; ii++) |
1008 cnt++; | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1009 |
10311 | 1010 ComplexColumnVector tmp (cnt); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1011 |
10311 | 1012 cnt = 0; |
1013 for (int ii = 0; ii < nn; ii++) | |
1014 tmp(cnt++) = xalpha(ii) / xbeta(ii); | |
1015 | |
1016 gev = tmp; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1017 } |
3185 | 1018 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1019 { |
3185 | 1020 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1021 std::cout << "qz: computing generalized eigenvalues" << std::endl; |
3185 | 1022 #endif |
3183 | 1023 |
10311 | 1024 // Return finite generalized eigenvalues. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1025 int cnt = 0; |
3185 | 1026 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1027 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1028 if (betar(ii) != 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1029 cnt++; |
3185 | 1030 |
10311 | 1031 ComplexColumnVector tmp (cnt); |
3185 | 1032 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1033 cnt = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1034 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1035 if (betar(ii) != 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1036 tmp(cnt++) = Complex(alphar(ii), alphai(ii))/betar(ii); |
10311 | 1037 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1038 gev = tmp; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1039 } |
3183 | 1040 } |
1041 | |
10311 | 1042 // Right, left eigenvector matrices. |
3185 | 1043 if (nargout >= 5) |
3183 | 1044 { |
10311 | 1045 // Which side to compute? |
1046 char side = (nargout == 5 ? 'R' : 'B'); | |
1047 // Compute all of them and backtransform | |
1048 char howmny = 'B'; | |
1049 // Dummy pointer; select is not used. | |
1050 octave_idx_type *select = 0; | |
3185 | 1051 |
1052 if (complex_case) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1053 { |
10311 | 1054 CVL = CQ; |
1055 CVR = CZ; | |
1056 ComplexRowVector cwork2 (2 * nn); | |
1057 RowVector rwork2 (8 * nn); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1058 octave_idx_type m; |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1059 |
10311 | 1060 F77_XFCN (ztgevc, ZTGEVC, |
1061 (F77_CONST_CHAR_ARG2 (&side, 1), | |
1062 F77_CONST_CHAR_ARG2 (&howmny, 1), | |
1063 select, nn, caa.fortran_vec (), nn, cbb.fortran_vec (), | |
1064 nn, CVL.fortran_vec (), nn, CVR.fortran_vec (), nn, nn, | |
1065 m, cwork2.fortran_vec (), rwork2.fortran_vec (), info | |
1066 F77_CHAR_ARG_LEN (1) | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
1067 F77_CHAR_ARG_LEN (1))); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1068 } |
3185 | 1069 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1070 { |
3185 | 1071 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1072 std::cout << "qz: computing generalized eigenvectors" << std::endl; |
3185 | 1073 #endif |
1074 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1075 VL = QQ; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1076 VR = ZZ; |
10311 | 1077 octave_idx_type m; |
1078 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1079 F77_XFCN (dtgevc, DTGEVC, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1080 (F77_CONST_CHAR_ARG2 (&side, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1081 F77_CONST_CHAR_ARG2 (&howmny, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1082 select, nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1083 nn, VL.fortran_vec (), nn, VR.fortran_vec (), nn, nn, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1084 m, work.fortran_vec (), info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1085 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1086 F77_CHAR_ARG_LEN (1))); |
3185 | 1087 |
10311 | 1088 // Now construct the complex form of VV, WW. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1089 int jj = 0; |
3185 | 1090 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1091 while (jj < nn) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1092 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1093 OCTAVE_QUIT; |
4153 | 1094 |
10311 | 1095 // See if real or complex eigenvalue. |
1096 | |
1097 // Column increment; assume complex eigenvalue. | |
1098 int cinc = 2; | |
3185 | 1099 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1100 if (jj == (nn-1)) |
10311 | 1101 // Single column. |
1102 cinc = 1; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1103 else if (aa(jj+1,jj) == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1104 cinc = 1; |
3185 | 1105 |
10311 | 1106 // Now copy the eigenvector (s) to CVR, CVL. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1107 if (cinc == 1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1108 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1109 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1110 CVR(ii,jj) = VR(ii,jj); |
3185 | 1111 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1112 if (side == 'B') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1113 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1114 CVL(ii,jj) = VL(ii,jj); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1115 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1116 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1117 { |
10311 | 1118 // Double column; complex vector. |
3185 | 1119 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1120 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1121 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1122 CVR(ii,jj) = Complex (VR(ii,jj), VR(ii,jj+1)); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1123 CVR(ii,jj+1) = Complex (VR(ii,jj), -VR(ii,jj+1)); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1124 } |
3183 | 1125 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1126 if (side == 'B') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1127 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1128 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1129 CVL(ii,jj) = Complex (VL(ii,jj), VL(ii,jj+1)); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1130 CVL(ii,jj+1) = Complex (VL(ii,jj), -VL(ii,jj+1)); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1131 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1132 } |
3185 | 1133 |
10311 | 1134 // Advance to next eigenvectors (if any). |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1135 jj += cinc; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1136 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1137 } |
3183 | 1138 } |
3185 | 1139 |
1140 switch (nargout) | |
1141 { | |
1142 case 7: | |
1143 retval(6) = gev; | |
1144 | |
10311 | 1145 case 6: |
1146 // Return eigenvectors. | |
3185 | 1147 retval(5) = CVL; |
1148 | |
10311 | 1149 case 5: |
1150 // Return eigenvectors. | |
3185 | 1151 retval(4) = CVR; |
1152 | |
1153 case 4: | |
1154 if (nargin == 3) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1155 { |
3185 | 1156 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1157 std::cout << "qz: sort: retval(3) = gev = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1158 octave_print_internal (std::cout, gev); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1159 std::cout << std::endl; |
3185 | 1160 #endif |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1161 retval(3) = gev; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1162 } |
3185 | 1163 else |
10311 | 1164 { |
1165 if (complex_case) | |
1166 retval(3) = CZ; | |
1167 else | |
1168 retval(3) = ZZ; | |
1169 } | |
3185 | 1170 |
1171 case 3: | |
1172 if (nargin == 3) | |
14876
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1173 { |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1174 if (complex_case) |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1175 retval(2) = CZ; |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1176 else |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1177 retval(2) = ZZ; |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1178 } |
3185 | 1179 else |
10311 | 1180 { |
1181 if (complex_case) | |
1182 retval(2) = CQ.hermitian (); | |
1183 else | |
1184 retval(2) = QQ.transpose (); | |
1185 } | |
1186 | |
3185 | 1187 case 2: |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1188 { |
10311 | 1189 if (complex_case) |
1190 { | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1191 #ifdef DEBUG |
14844
5bc9b9cb4362
maint: Use Octave coding conventions for cuddled parenthesis in retval assignments.
Rik <octave@nomad.inbox5.com>
parents:
14501
diff
changeset
|
1192 std::cout << "qz: retval(1) = cbb = " << std::endl; |
10311 | 1193 octave_print_internal (std::cout, cbb, 0); |
1194 std::cout << std::endl << "qz: retval(0) = caa = " <<std::endl; | |
1195 octave_print_internal (std::cout, caa, 0); | |
1196 std::cout << std::endl; | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
1197 #endif |
10311 | 1198 retval(1) = cbb; |
1199 retval(0) = caa; | |
1200 } | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1201 else |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1202 { |
3185 | 1203 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1204 std::cout << "qz: retval(1) = bb = " << std::endl; |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1205 octave_print_internal (std::cout, bb, 0); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1206 std::cout << std::endl << "qz: retval(0) = aa = " <<std::endl; |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1207 octave_print_internal (std::cout, aa, 0); |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1208 std::cout << std::endl; |
3185 | 1209 #endif |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1210 retval(1) = bb; |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1211 retval(0) = aa; |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1212 } |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
1213 } |
10311 | 1214 break; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1215 |
3185 | 1216 |
1217 case 1: | |
1218 case 0: | |
1219 #ifdef DEBUG | |
3538 | 1220 std::cout << "qz: retval(0) = gev = " << gev << std::endl; |
3185 | 1221 #endif |
1222 retval(0) = gev; | |
1223 break; | |
1224 | |
1225 default: | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
1226 error ("qz: too many return arguments"); |
3185 | 1227 break; |
3183 | 1228 } |
1229 | |
3185 | 1230 #ifdef DEBUG |
3538 | 1231 std::cout << "qz: exiting (at long last)" << std::endl; |
3185 | 1232 #endif |
3183 | 1233 |
1234 return retval; | |
1235 } | |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1236 |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1237 /* |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1238 %!shared a, b, c |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1239 %! a = [1 2; 0 3]; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1240 %! b = [1 0; 0 0]; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1241 %! c = [0 1; 0 0]; |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1242 %!assert (qz (a,b), 1) |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1243 %!assert (isempty (qz (a,c))) |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1244 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1245 ## Exaple 7.7.3 in Golub & Van Loan |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1246 %!test |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1247 %! a = [ 10 1 2; |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1248 %! 1 2 -1; |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1249 %! 1 1 2]; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1250 %! b = reshape (1:9,3,3); |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1251 %! [aa, bb, q, z, v, w, lambda] = qz (a, b); |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1252 %! sz = length (lambda); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1253 %! observed = (b * v * diag ([lambda;0])) (:, 1:sz); |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
1254 %! assert ((a*v)(:, 1:sz), observed, norm (observed) * 1e-14); |
13088
4ffea87ad71b
codesprint: Fix tolerance for qz.cc tests.
Ben Abbott <bpabbott@mac.com>
parents:
13074
diff
changeset
|
1255 %! observed = (diag ([lambda;0]) * w' * b) (1:sz, :); |
18678
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
1256 %! assert ((w'*a)(1:sz, :) , observed, norm (observed) * 1e-13); |
13088
4ffea87ad71b
codesprint: Fix tolerance for qz.cc tests.
Ben Abbott <bpabbott@mac.com>
parents:
13074
diff
changeset
|
1257 %! 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
|
1258 %! assert (q * b * z, bb, norm (bb) * 1e-14); |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1259 |
14876
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1260 %!test |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1261 %! A = [0, 0, -1, 0; 1, 0, 0, 0; -1, 0, -2, -1; 0, -1, 1, 0]; |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1262 %! B = [0, 0, 0, 0; 0, 1, 0, 0; 0, 0, 1, 0; 0, 0, 0, 1]; |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1263 %! [AA, BB, Q, Z1] = qz (A, B); |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1264 %! [AA, BB, Z2] = qz (A, B, '-'); |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1265 %! assert (Z1, Z2); |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1266 */ |