Mercurial > octave-nkf
annotate src/DLD-FUNCTIONS/qz.cc @ 10791:3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
interpreter/doccheck: New directory for spelling/grammar scripts.
interpreter/doccheck/README: Instructions for using scripts.
interpreter/doccheck/spellcheck: Script to spellcheck a Texinfo file.
interpreter/doccheck/aspell.conf: GNU Aspell configuration file for
Octave documentation.
interpreter/doccheck/aspell-octave.en.pws: Private Aspell dictionary.
interpreter/doccheck/add_to_aspell_dict: Script to add new
Octave-specific words to
private Aspell dictionary.
interpreter/octave.texi: New @nospell macro which forces Aspell
to ignore the word marked by the macro.
interpreter/mk_doc_cache.m: Skip new @nospell macro when building
doc_cache.
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Sat, 17 Jul 2010 19:53:01 -0700 |
parents | 12884915a8e4 |
children | 89f4d7e294cc |
rev | line source |
---|---|
3183 | 1 /* |
2 | |
8920 | 3 Copyright (C) 1998, 1999, 2000, 2002, 2003, 2004, 2005, 2006, 2007, |
4 2008, 2009 A. S. Hodel | |
3183 | 5 |
6 This file is part of Octave. | |
7 | |
8 Octave is free software; you can redistribute it and/or modify it | |
9 under the terms of the GNU General Public License as published by the | |
7016 | 10 Free Software Foundation; either version 3 of the License, or (at your |
11 option) any later version. | |
3183 | 12 |
13 Octave is distributed in the hope that it will be useful, but WITHOUT | |
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
16 for more details. | |
17 | |
18 You should have received a copy of the GNU General Public License | |
7016 | 19 along with Octave; see the file COPYING. If not, see |
20 <http://www.gnu.org/licenses/>. | |
3183 | 21 |
22 */ | |
23 | |
24 // Generalized eigenvalue balancing via LAPACK | |
3911 | 25 |
26 // Author: A. S. Hodel <scotte@eng.auburn.edu> | |
3183 | 27 |
28 #undef DEBUG | |
29 #undef DEBUG_SORT | |
30 #undef DEBUG_EIG | |
31 | |
9786
2c279308f6ab
fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
32 #ifdef HAVE_CONFIG_H |
2c279308f6ab
fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
33 #include <config.h> |
2c279308f6ab
fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents:
9758
diff
changeset
|
34 #endif |
3183 | 35 |
36 #include <cfloat> | |
4051 | 37 |
3523 | 38 #include <iostream> |
4051 | 39 #include <iomanip> |
3183 | 40 |
41 #include "CmplxQRP.h" | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
42 #include "CmplxQR.h" |
3183 | 43 #include "dbleQR.h" |
4153 | 44 #include "f77-fcn.h" |
7231 | 45 #include "lo-math.h" |
4153 | 46 #include "quit.h" |
47 | |
3183 | 48 #include "defun-dld.h" |
49 #include "error.h" | |
50 #include "gripes.h" | |
51 #include "oct-obj.h" | |
52 #include "oct-map.h" | |
53 #include "ov.h" | |
54 #include "pager.h" | |
3185 | 55 #if defined (DEBUG) || defined (DEBUG_SORT) |
3183 | 56 #include "pr-output.h" |
57 #endif | |
58 #include "symtab.h" | |
59 #include "utils.h" | |
60 #include "variables.h" | |
61 | |
5275 | 62 typedef octave_idx_type (*sort_function) (const octave_idx_type& LSIZE, const double& ALPHA, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
63 const double& BETA, const double& S, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
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 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
78 F77_RET_T |
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 | |
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 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
100 F77_RET_T |
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 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
125 F77_RET_T |
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 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
156 F77_RET_T |
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, | |
166 const octave_idx_type& LDQ, | |
167 Complex* CZ, const octave_idx_type& LDZ, | |
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 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
207 F77_RET_T |
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 | |
291 DEFUN_DLD (qz, args, nargout, | |
3372 | 292 "-*- texinfo -*-\n\ |
293 @deftypefn {Loadable Function} {@var{lambda} =} qz (@var{a}, @var{b})\n\ | |
294 Generalized eigenvalue problem @math{A x = s B x},\n\ | |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10350
diff
changeset
|
295 QZ decomposition. There are three ways to call this function:\n\ |
3372 | 296 @enumerate\n\ |
297 @item @code{lambda = qz(A,B)}\n\ | |
3185 | 298 \n\ |
5016 | 299 Computes the generalized eigenvalues\n\ |
300 @tex\n\ | |
301 $\\lambda$\n\ | |
302 @end tex\n\ | |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
303 @ifnottex\n\ |
5016 | 304 @var{lambda}\n\ |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
305 @end ifnottex\n\ |
5016 | 306 of @math{(A - s B)}.\n\ |
3500 | 307 @item @code{[AA, BB, Q, Z, V, W, lambda] = qz (A, B)}\n\ |
3185 | 308 \n\ |
3372 | 309 Computes qz decomposition, generalized eigenvectors, and \n\ |
5481 | 310 generalized eigenvalues of @math{(A - sB)}\n\ |
5016 | 311 @tex\n\ |
312 $$ AV = BV{ \\rm diag }(\\lambda) $$\n\ | |
313 $$ W^T A = { \\rm diag }(\\lambda)W^T B $$\n\ | |
314 $$ AA = Q^T AZ, BB = Q^T BZ $$\n\ | |
315 @end tex\n\ | |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
316 @ifnottex\n\ |
3372 | 317 @example\n\ |
318 @group\n\ | |
5481 | 319 \n\ |
10311 | 320 A * V = B * V * diag (lambda)\n\ |
321 W' * A = diag (lambda) * W' * B\n\ | |
322 AA = Q * A * Z, BB = Q * B * Z\n\ | |
5481 | 323 \n\ |
3372 | 324 @end group\n\ |
325 @end example\n\ | |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
326 @end ifnottex\n\ |
5016 | 327 with @var{Q} and @var{Z} orthogonal (unitary)= @var{I}\n\ |
3185 | 328 \n\ |
5016 | 329 @item @code{[AA,BB,Z@{, lambda@}] = qz(A,B,opt)}\n\ |
3185 | 330 \n\ |
3372 | 331 As in form [2], but allows ordering of generalized eigenpairs\n\ |
5481 | 332 for (e.g.) solution of discrete time algebraic Riccati equations.\n\ |
333 Form 3 is not available for complex matrices, and does not compute\n\ | |
334 the generalized 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
|
335 \n\ |
3372 | 336 @table @var\n\ |
337 @item opt\n\ | |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8927
diff
changeset
|
338 for ordering eigenvalues of the GEP pencil. The leading block\n\ |
5481 | 339 of the revised pencil contains all eigenvalues that satisfy:\n\ |
3372 | 340 @table @code\n\ |
341 @item \"N\"\n\ | |
5481 | 342 = unordered (default) \n\ |
3183 | 343 \n\ |
3372 | 344 @item \"S\"\n\ |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10350
diff
changeset
|
345 = small: leading block has all |lambda| <= 1 \n\ |
3185 | 346 \n\ |
3372 | 347 @item \"B\"\n\ |
5481 | 348 = big: leading block has all |lambda| >= 1 \n\ |
3372 | 349 \n\ |
350 @item \"-\"\n\ | |
5481 | 351 = negative real part: leading block has all eigenvalues\n\ |
352 in the open left half-plane\n\ | |
3372 | 353 \n\ |
354 @item \"+\"\n\ | |
9209
923c7cb7f13f
Simplify TeXinfo files by eliminating redundant @iftex followed by @tex construction.
Rik <rdrider0-list@yahoo.com>
parents:
9064
diff
changeset
|
355 = non-negative real part: leading block has all eigenvalues\n\ |
5481 | 356 in the closed right half-plane\n\ |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8927
diff
changeset
|
357 @end table\n\ |
3372 | 358 @end table\n\ |
359 @end enumerate\n\ | |
3183 | 360 \n\ |
3372 | 361 Note: qz performs permutation balancing, but not scaling (see balance).\n\ |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10350
diff
changeset
|
362 The order of output arguments was selected for compatibility with @sc{matlab}\n\ |
3183 | 363 \n\ |
8286
6f2d95255911
fix @seealso references to point to existing anchors
Thorsten Meyer <thorsten.meyier@gmx.de>
parents:
7520
diff
changeset
|
364 @seealso{balance, eig, schur}\n\ |
3372 | 365 @end deftypefn") |
3183 | 366 { |
367 octave_value_list retval; | |
368 int nargin = args.length (); | |
369 | |
3185 | 370 #ifdef DEBUG |
3538 | 371 std::cout << "qz: nargin = " << nargin << ", nargout = " << nargout << std::endl; |
3185 | 372 #endif |
3183 | 373 |
3185 | 374 if (nargin < 2 || nargin > 3 || nargout > 7) |
375 { | |
5823 | 376 print_usage (); |
3185 | 377 return retval; |
378 } | |
379 else if (nargin == 3 && (nargout < 3 || nargout > 4)) | |
380 { | |
3427 | 381 error ("qz: invalid number of output arguments for form [3] call"); |
3185 | 382 return retval; |
383 } | |
3183 | 384 |
3185 | 385 #ifdef DEBUG |
3538 | 386 std::cout << "qz: determine ordering option" << std::endl; |
3185 | 387 #endif |
3183 | 388 |
10311 | 389 // 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
|
390 volatile char ord_job = 0; |
3183 | 391 static double safmin; |
3185 | 392 |
393 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
|
394 ord_job = 'N'; |
3185 | 395 else if (!args(2).is_string ()) |
396 { | |
397 error ("qz: argument 3 must be a string"); | |
398 return retval; | |
399 } | |
400 else | |
401 { | |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
402 std::string tmp = args(2).string_value (); |
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
403 |
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
404 if (! tmp.empty ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
405 ord_job = tmp[0]; |
3183 | 406 |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
407 if (! (ord_job == 'N' || ord_job == 'n' |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
408 || ord_job == 'S' || ord_job == 's' |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
409 || ord_job == 'B' || ord_job == 'b' |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
410 || ord_job == '+' || ord_job == '-')) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
411 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
412 error ("qz: invalid order option"); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
413 return retval; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
414 } |
3185 | 415 |
416 // overflow constant required by dlag2 | |
4552 | 417 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
|
418 safmin |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
419 F77_CHAR_ARG_LEN (1)); |
3183 | 420 |
3185 | 421 #ifdef DEBUG_EIG |
3538 | 422 std::cout << "qz: initial value of safmin=" << setiosflags (std::ios::scientific) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
423 << safmin << std::endl; |
3185 | 424 #endif |
3183 | 425 |
10311 | 426 // Some machines (e.g., DEC alpha) get safmin = 0; |
427 // for these, use eps instead to avoid problems in dlag2. | |
3185 | 428 if (safmin == 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
429 { |
3185 | 430 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
431 std::cout << "qz: DANGER WILL ROBINSON: safmin is 0!" << std::endl; |
3185 | 432 #endif |
3183 | 433 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
434 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
|
435 safmin |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
436 F77_CHAR_ARG_LEN (1)); |
3185 | 437 |
438 #ifdef DEBUG_EIG | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
439 std::cout << "qz: safmin set to " << setiosflags (std::ios::scientific) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
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 | |
3531 | 453 std::cout << "argument 1 dimensions: (" << nn << "," << args(0).columns () << ")" |
3538 | 454 << std::endl; |
3185 | 455 #endif |
456 | |
457 int arg_is_empty = empty_arg ("qz", nn, args(0).columns ()); | |
458 | |
3183 | 459 if (arg_is_empty < 0) |
3185 | 460 { |
461 gripe_empty_arg ("qz: parameter 1", 0); | |
462 return retval; | |
463 } | |
3183 | 464 else if (arg_is_empty > 0) |
3185 | 465 { |
466 gripe_empty_arg ("qz: parameter 1; continuing", 0); | |
467 return octave_value_list (2, Matrix ()); | |
468 } | |
469 else if (args(0).columns () != nn) | |
470 { | |
471 gripe_square_matrix_required ("qz"); | |
472 return retval; | |
473 } | |
3183 | 474 |
10311 | 475 // Argument 1: dimensions look good; get the value. |
3183 | 476 Matrix aa; |
477 ComplexMatrix caa; | |
3185 | 478 |
479 if (args(0).is_complex_type ()) | |
3183 | 480 caa = args(0).complex_matrix_value (); |
3185 | 481 else |
3183 | 482 aa = args(0).matrix_value (); |
3185 | 483 |
484 if (error_state) | |
3183 | 485 return retval; |
486 | |
3185 | 487 #ifdef DEBUG |
3538 | 488 std::cout << "qz: check argument 2" << std::endl; |
3185 | 489 #endif |
3183 | 490 |
10311 | 491 // Extract argument 2 (bb, or cbb if complex). |
3185 | 492 if ((nn != args(1).columns ()) || (nn != args(1).rows ())) |
493 { | |
494 gripe_nonconformant (); | |
495 return retval; | |
496 } | |
497 | |
3183 | 498 Matrix bb; |
499 ComplexMatrix cbb; | |
3185 | 500 |
501 if (args(1).is_complex_type ()) | |
3183 | 502 cbb = args(1).complex_matrix_value (); |
503 else | |
504 bb = args(1).matrix_value (); | |
3185 | 505 |
506 if (error_state) | |
3183 | 507 return retval; |
508 | |
509 // Both matrices loaded, now let's check what kind of arithmetic: | |
10311 | 510 // declared volatile to avoid compiler warnings about long jumps, |
511 // vforks. | |
3185 | 512 |
10311 | 513 volatile int complex_case |
3185 | 514 = (args(0).is_complex_type () || args(1).is_complex_type ()); |
10311 | 515 |
3185 | 516 if (nargin == 3 && complex_case) |
517 { | |
518 error ("qz: cannot re-order complex qz decomposition."); | |
519 return retval; | |
520 } | |
3183 | 521 |
10311 | 522 // First, declare variables used in both the real and complex case. |
3183 | 523 Matrix QQ(nn,nn), ZZ(nn,nn), VR(nn,nn), VL(nn,nn); |
524 RowVector alphar(nn), alphai(nn), betar(nn); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
525 ComplexRowVector xalpha(nn), xbeta(nn); |
3185 | 526 ComplexMatrix CQ(nn,nn), CZ(nn,nn), CVR(nn,nn), CVL(nn,nn); |
5275 | 527 octave_idx_type ilo, ihi, info; |
3185 | 528 char compq = (nargout >= 3 ? 'V' : 'N'); |
529 char compz = (nargout >= 4 ? 'V' : 'N'); | |
3183 | 530 |
10311 | 531 // Initialize Q, Z to identity if we need either of them. |
3185 | 532 if (compq == 'V' || compz == 'V') |
5275 | 533 for (octave_idx_type ii = 0; ii < nn; ii++) |
534 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
|
535 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
536 OCTAVE_QUIT; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
537 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
|
538 } |
3183 | 539 |
10311 | 540 // Always perform permutation balancing. |
4552 | 541 const char bal_job = 'P'; |
10311 | 542 RowVector lscale (nn), rscale (nn), work (6 * nn), rwork (nn); |
3183 | 543 |
3185 | 544 if (complex_case) |
545 { | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
546 #ifdef DEBUG |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
547 if (compq == 'V') |
10311 | 548 std::cout << "qz: performing balancing; CQ=" << std::endl << CQ << std::endl; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
549 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
550 if (args(0).is_real_type ()) |
10311 | 551 caa = ComplexMatrix (aa); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
552 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
553 if (args(1).is_real_type ()) |
10311 | 554 cbb = ComplexMatrix (bb); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
555 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
556 if (compq == 'V') |
10311 | 557 CQ = ComplexMatrix (QQ); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
558 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
559 if (compz == 'V') |
10311 | 560 CZ = ComplexMatrix (ZZ); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
561 |
10311 | 562 F77_XFCN (zggbal, ZGGBAL, |
563 (F77_CONST_CHAR_ARG2 (&bal_job, 1), | |
564 nn, caa.fortran_vec (), nn, cbb.fortran_vec (), | |
565 nn, ilo, ihi, lscale.fortran_vec (), | |
566 rscale.fortran_vec (), work.fortran_vec (), info | |
567 F77_CHAR_ARG_LEN (1))); | |
3185 | 568 } |
3183 | 569 else |
3185 | 570 { |
571 #ifdef DEBUG | |
572 if (compq == 'V') | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
573 std::cout << "qz: performing balancing; QQ=" << std::endl << QQ << std::endl; |
3185 | 574 #endif |
3183 | 575 |
3185 | 576 F77_XFCN (dggbal, DGGBAL, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
577 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
578 nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
579 nn, ilo, ihi, lscale.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
580 rscale.fortran_vec (), work.fortran_vec (), info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
581 F77_CHAR_ARG_LEN (1))); |
3185 | 582 } |
3183 | 583 |
584 // Since we just want the balancing matrices, we can use dggbal | |
10311 | 585 // for both the real and complex cases; left first |
3185 | 586 |
10311 | 587 #if 0 |
588 if (compq == 'V') | |
3185 | 589 { |
590 F77_XFCN (dggbak, DGGBAK, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
591 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
592 F77_CONST_CHAR_ARG2 ("L", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
593 nn, ilo, ihi, lscale.data (), rscale.data (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
594 nn, QQ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
595 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
596 F77_CHAR_ARG_LEN (1))); |
3183 | 597 |
3185 | 598 #ifdef DEBUG |
599 if (compq == 'V') | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
600 std::cout << "qz: balancing done; QQ=" << std::endl << QQ << std::endl; |
3185 | 601 #endif |
3183 | 602 } |
603 | |
604 // then right | |
3185 | 605 if (compz == 'V') |
606 { | |
4552 | 607 F77_XFCN (dggbak, DGGBAK, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
608 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
609 F77_CONST_CHAR_ARG2 ("R", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
610 nn, ilo, ihi, lscale.data (), rscale.data (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
611 nn, ZZ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
612 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
613 F77_CHAR_ARG_LEN (1))); |
3183 | 614 |
3185 | 615 #ifdef DEBUG |
616 if (compz == 'V') | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
617 std::cout << "qz: balancing done; ZZ=" << std::endl << ZZ << std::endl; |
3185 | 618 #endif |
10311 | 619 } |
620 #endif | |
3183 | 621 |
622 static char qz_job; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
623 qz_job = (nargout < 2 ? 'E' : 'S'); |
3185 | 624 |
3183 | 625 if (complex_case) |
3185 | 626 { |
10311 | 627 // Complex case. |
628 | |
629 // The QR decomposition of cbb. | |
630 ComplexQR cbqr (cbb); | |
631 // The R matrix of QR decomposition for cbb. | |
632 cbb = cbqr.R (); | |
633 // (Q*)caa for following work. | |
634 caa = (cbqr.Q ().hermitian ()) * caa; | |
635 CQ = CQ * cbqr.Q (); | |
636 | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
637 F77_XFCN (zgghrd, ZGGHRD, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
638 (F77_CONST_CHAR_ARG2 (&compq, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
639 F77_CONST_CHAR_ARG2 (&compz, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
640 nn, ilo, ihi, caa.fortran_vec (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
641 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
|
642 CZ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
643 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
644 F77_CHAR_ARG_LEN (1))); |
10311 | 645 |
646 ComplexRowVector cwork (1 * nn); | |
647 | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
648 F77_XFCN (zhgeqz, ZHGEQZ, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
649 (F77_CONST_CHAR_ARG2 (&qz_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
650 F77_CONST_CHAR_ARG2 (&compq, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
651 F77_CONST_CHAR_ARG2 (&compz, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
652 nn, ilo, ihi, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
653 caa.fortran_vec (), nn, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
654 cbb.fortran_vec (),nn, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
655 xalpha.fortran_vec (), xbeta.fortran_vec (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
656 CQ.fortran_vec (), nn, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
657 CZ.fortran_vec (), nn, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
658 cwork.fortran_vec (), nn, rwork.fortran_vec (), info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
659 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
660 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
661 F77_CHAR_ARG_LEN (1))); |
3185 | 662 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
663 if (compq == 'V') |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
664 { |
10311 | 665 // Left eigenvector. |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
666 F77_XFCN (zggbak, ZGGBAK, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
667 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
668 F77_CONST_CHAR_ARG2 ("L", 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
669 nn, ilo, ihi, lscale.data (), rscale.data (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
670 nn, CQ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
671 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
672 F77_CHAR_ARG_LEN (1))); |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
673 } |
3185 | 674 |
10311 | 675 // Right eigenvector. |
3185 | 676 if (compz == 'V') |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
677 { |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
678 F77_XFCN (zggbak, ZGGBAK, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
679 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
680 F77_CONST_CHAR_ARG2 ("R", 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
681 nn, ilo, ihi, lscale.data (), rscale.data (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
682 nn, CZ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
683 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
684 F77_CHAR_ARG_LEN (1))); |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
685 } |
3185 | 686 |
687 } | |
10311 | 688 else |
3185 | 689 { |
690 #ifdef DEBUG | |
3538 | 691 std::cout << "qz: peforming qr decomposition of bb" << std::endl; |
3185 | 692 #endif |
3183 | 693 |
10311 | 694 // Compute the QR factorization of bb. |
3185 | 695 QR bqr (bb); |
696 | |
697 #ifdef DEBUG | |
3538 | 698 std::cout << "qz: qr (bb) done; now peforming qz decomposition" << std::endl; |
3185 | 699 #endif |
3183 | 700 |
3185 | 701 bb = bqr.R (); |
702 | |
703 #ifdef DEBUG | |
3538 | 704 std::cout << "qz: extracted bb" << std::endl; |
3185 | 705 #endif |
3183 | 706 |
10311 | 707 aa = (bqr.Q ()).transpose () * aa; |
3185 | 708 |
709 #ifdef DEBUG | |
3538 | 710 std::cout << "qz: updated aa " << std::endl; |
711 std::cout << "bqr.Q () = " << std::endl << bqr.Q () << std::endl; | |
3183 | 712 |
3185 | 713 if (compq == 'V') |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
714 std::cout << "QQ =" << QQ << std::endl; |
3185 | 715 #endif |
3183 | 716 |
3185 | 717 if (compq == 'V') |
10311 | 718 QQ = QQ * bqr.Q (); |
3183 | 719 |
3185 | 720 #ifdef DEBUG |
3538 | 721 std::cout << "qz: precursors done..." << std::endl; |
3185 | 722 #endif |
3183 | 723 |
3185 | 724 #ifdef DEBUG |
3538 | 725 std::cout << "qz: compq = " << compq << ", compz = " << compz << std::endl; |
3185 | 726 #endif |
3183 | 727 |
10311 | 728 // Reduce to generalized hessenberg form. |
3185 | 729 F77_XFCN (dgghrd, DGGHRD, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
730 (F77_CONST_CHAR_ARG2 (&compq, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
731 F77_CONST_CHAR_ARG2 (&compz, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
732 nn, ilo, ihi, aa.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
733 nn, bb.fortran_vec (), nn, QQ.fortran_vec (), nn, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
734 ZZ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
735 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
736 F77_CHAR_ARG_LEN (1))); |
3183 | 737 |
10311 | 738 // Check if just computing generalized eigenvalues or if we're |
739 // actually computing the decomposition. | |
3183 | 740 |
10311 | 741 // Reduce to generalized Schur form. |
3185 | 742 F77_XFCN (dhgeqz, DHGEQZ, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
743 (F77_CONST_CHAR_ARG2 (&qz_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
744 F77_CONST_CHAR_ARG2 (&compq, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
745 F77_CONST_CHAR_ARG2 (&compz, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
746 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
|
747 nn, alphar.fortran_vec (), alphai.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
748 betar.fortran_vec (), QQ.fortran_vec (), nn, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
749 ZZ.fortran_vec (), nn, work.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
750 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
751 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
752 F77_CHAR_ARG_LEN (1))); |
10311 | 753 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
754 if (compq == 'V') |
10311 | 755 { |
756 F77_XFCN (dggbak, DGGBAK, | |
757 (F77_CONST_CHAR_ARG2 (&bal_job, 1), | |
758 F77_CONST_CHAR_ARG2 ("L", 1), | |
759 nn, ilo, ihi, lscale.data (), rscale.data (), | |
760 nn, QQ.fortran_vec (), nn, info | |
761 F77_CHAR_ARG_LEN (1) | |
762 F77_CHAR_ARG_LEN (1))); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
763 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
764 #ifdef DEBUG |
10311 | 765 if (compq == 'V') |
766 std::cout << "qz: balancing done; QQ=" << std::endl << QQ << std::endl; | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
767 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
768 } |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
769 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
770 // then right |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
771 if (compz == 'V') |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
772 { |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
773 F77_XFCN (dggbak, DGGBAK, |
10311 | 774 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
775 F77_CONST_CHAR_ARG2 ("R", 1), | |
776 nn, ilo, ihi, lscale.data (), rscale.data (), | |
777 nn, ZZ.fortran_vec (), nn, info | |
778 F77_CHAR_ARG_LEN (1) | |
779 F77_CHAR_ARG_LEN (1))); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
780 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
781 #ifdef DEBUG |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
782 if (compz == 'V') |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
783 std::cout << "qz: balancing done; ZZ=" << std::endl << ZZ << std::endl; |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
784 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
785 } |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
786 |
3185 | 787 } |
3183 | 788 |
10311 | 789 // 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
|
790 if (! (ord_job == 'N' || ord_job == 'n')) |
3183 | 791 { |
3185 | 792 if (complex_case) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
793 { |
10311 | 794 // Probably not needed, but better be safe. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
795 error ("qz: cannot re-order complex qz decomposition."); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
796 return retval; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
797 } |
3185 | 798 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
799 { |
3185 | 800 #ifdef DEBUG_SORT |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
801 std::cout << "qz: ordering eigenvalues: ord_job = " |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
802 << ord_job << std::endl; |
3185 | 803 #endif |
3183 | 804 |
10311 | 805 // Declared static to avoid vfork/long jump compiler complaints. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
806 static sort_function sort_test; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
807 sort_test = 0; |
3183 | 808 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
809 switch (ord_job) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
810 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
811 case 'S': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
812 case 's': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
813 sort_test = &fin; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
814 break; |
3183 | 815 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
816 case 'B': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
817 case 'b': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
818 sort_test = &fout; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
819 break; |
3183 | 820 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
821 case '+': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
822 sort_test = &fcrhp; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
823 break; |
3183 | 824 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
825 case '-': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
826 sort_test = &folhp; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
827 break; |
3185 | 828 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
829 default: |
10311 | 830 // Invalid order option (should never happen, since we |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
831 // checked the options at the top). |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
832 panic_impossible (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
833 break; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
834 } |
3183 | 835 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
836 octave_idx_type ndim, fail; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
837 double inf_norm; |
3185 | 838 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
839 F77_XFCN (xdlange, XDLANGE, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
840 (F77_CONST_CHAR_ARG2 ("I", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
841 nn, nn, aa.data (), nn, work.fortran_vec (), inf_norm |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
842 F77_CHAR_ARG_LEN (1))); |
3185 | 843 |
10311 | 844 double eps = DBL_EPSILON * inf_norm * nn; |
3185 | 845 |
846 #ifdef DEBUG_SORT | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
847 std::cout << "qz: calling dsubsp: aa=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
848 octave_print_internal (std::cout, aa, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
849 std::cout << std::endl << "bb=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
850 octave_print_internal (std::cout, bb, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
851 if (compz == 'V') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
852 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
853 std::cout << std::endl << "ZZ=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
854 octave_print_internal (std::cout, ZZ, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
855 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
856 std::cout << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
857 std::cout << "alphar = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
858 octave_print_internal (std::cout, (Matrix) alphar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
859 std::cout << std::endl << "alphai = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
860 octave_print_internal (std::cout, (Matrix) alphai, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
861 std::cout << std::endl << "beta = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
862 octave_print_internal (std::cout, (Matrix) betar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
863 std::cout << std::endl; |
3185 | 864 #endif |
865 | |
10350
12884915a8e4
merge MArray classes & improve Array interface
Jaroslav Hajek <highegg@gmail.com>
parents:
10311
diff
changeset
|
866 Array<octave_idx_type> ind (nn, 1); |
3550 | 867 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
868 F77_XFCN (dsubsp, DSUBSP, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
869 (nn, nn, aa.fortran_vec (), bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
870 ZZ.fortran_vec (), sort_test, eps, ndim, fail, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
871 ind.fortran_vec ())); |
3185 | 872 |
873 #ifdef DEBUG | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
874 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
875 octave_print_internal (std::cout, aa, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
876 std::cout << std::endl << "bb=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
877 octave_print_internal (std::cout, bb, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
878 if (compz == 'V') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
879 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
880 std::cout << std::endl << "ZZ=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
881 octave_print_internal (std::cout, ZZ, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
882 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
883 std::cout << std::endl; |
3185 | 884 #endif |
885 | |
10311 | 886 // Manually update alphar, alphai, betar. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
887 static int jj; |
3185 | 888 |
10311 | 889 jj = 0; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
890 while (jj < nn) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
891 { |
3185 | 892 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
893 std::cout << "computing gen eig #" << jj << std::endl; |
3185 | 894 #endif |
895 | |
10311 | 896 // Number of zeros in this block. |
897 static int zcnt; | |
3185 | 898 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
899 if (jj == (nn-1)) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
900 zcnt = 1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
901 else if (aa(jj+1,jj) == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
902 zcnt = 1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
903 else zcnt = 2; |
3185 | 904 |
10311 | 905 if (zcnt == 1) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
906 { |
10311 | 907 // Real zero. |
3185 | 908 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
909 std::cout << " single gen eig:" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
910 std::cout << " alphar(" << jj << ") = " << aa(jj,jj) << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
911 std::cout << " betar( " << jj << ") = " << bb(jj,jj) << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
912 std::cout << " alphai(" << jj << ") = 0" << std::endl; |
3185 | 913 #endif |
914 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
915 alphar(jj) = aa(jj,jj); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
916 alphai(jj) = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
917 betar(jj) = bb(jj,jj); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
918 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
919 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
920 { |
10311 | 921 // Complex conjugate pair. |
3185 | 922 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
923 std::cout << "qz: calling dlag2:" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
924 std::cout << "safmin=" |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
925 << setiosflags (std::ios::scientific) << safmin << std::endl; |
3185 | 926 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
927 for (int idr = jj; idr <= jj+1; idr++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
928 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
929 for (int idc = jj; idc <= jj+1; idc++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
930 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
931 std::cout << "aa(" << idr << "," << idc << ")=" |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
932 << aa(idr,idc) << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
933 std::cout << "bb(" << idr << "," << idc << ")=" |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
934 << bb(idr,idc) << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
935 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
936 } |
3185 | 937 #endif |
938 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
939 // FIXME -- probably should be using |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
940 // fortran_vec instead of &aa(jj,jj) here. |
4566 | 941 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
942 double scale1, scale2, wr1, wr2, wi; |
10311 | 943 const double *aa_ptr = aa.data () + jj * nn + jj; |
944 const double *bb_ptr = bb.data () + jj * nn + jj; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
945 F77_XFCN (dlag2, DLAG2, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
946 (aa_ptr, nn, bb_ptr, nn, safmin, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
947 scale1, scale2, wr1, wr2, wi)); |
3185 | 948 |
949 #ifdef DEBUG_EIG | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
950 std::cout << "dlag2 returns: scale1=" << scale1 |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
951 << "\tscale2=" << scale2 << std::endl |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
952 << "\twr1=" << wr1 << "\twr2=" << wr2 |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
953 << "\twi=" << wi << std::endl; |
3185 | 954 #endif |
955 | |
10311 | 956 // Just to be safe, check if it's a real pair. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
957 if (wi == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
958 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
959 alphar(jj) = wr1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
960 alphai(jj) = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
961 betar(jj) = scale1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
962 alphar(jj+1) = wr2; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
963 alphai(jj+1) = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
964 betar(jj+1) = scale2; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
965 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
966 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
967 { |
10311 | 968 alphar(jj) = alphar(jj+1) = wr1; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
969 alphai(jj) = -(alphai(jj+1) = wi); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
970 betar(jj) = betar(jj+1) = scale1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
971 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
972 } |
3185 | 973 |
10311 | 974 // Advance past this block. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
975 jj += zcnt; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
976 } |
3185 | 977 |
978 #ifdef DEBUG_SORT | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
979 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
980 octave_print_internal (std::cout, aa, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
981 std::cout << std::endl << "bb=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
982 octave_print_internal (std::cout, bb, 0); |
3185 | 983 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
984 if (compz == 'V') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
985 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
986 std::cout << std::endl << "ZZ=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
987 octave_print_internal (std::cout, ZZ, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
988 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
989 std::cout << std::endl << "qz: ndim=" << ndim << std::endl |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
990 << "fail=" << fail << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
991 std::cout << "alphar = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
992 octave_print_internal (std::cout, (Matrix) alphar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
993 std::cout << std::endl << "alphai = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
994 octave_print_internal (std::cout, (Matrix) alphai, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
995 std::cout << std::endl << "beta = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
996 octave_print_internal (std::cout, (Matrix) betar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
997 std::cout << std::endl; |
3185 | 998 #endif |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
999 } |
3183 | 1000 } |
3185 | 1001 |
10311 | 1002 // Compute generalized eigenvalues? |
3183 | 1003 ComplexColumnVector gev; |
3185 | 1004 |
1005 if (nargout < 2 || nargout == 7 || (nargin == 3 && nargout == 4)) | |
3183 | 1006 { |
3185 | 1007 if (complex_case) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1008 { |
10311 | 1009 int cnt = 0; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1010 |
10311 | 1011 for (int ii = 0; ii < nn; ii++) |
1012 cnt++; | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1013 |
10311 | 1014 ComplexColumnVector tmp (cnt); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1015 |
10311 | 1016 cnt = 0; |
1017 for (int ii = 0; ii < nn; ii++) | |
1018 tmp(cnt++) = xalpha(ii) / xbeta(ii); | |
1019 | |
1020 gev = tmp; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1021 } |
3185 | 1022 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1023 { |
3185 | 1024 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1025 std::cout << "qz: computing generalized eigenvalues" << std::endl; |
3185 | 1026 #endif |
3183 | 1027 |
10311 | 1028 // Return finite generalized eigenvalues. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1029 int cnt = 0; |
3185 | 1030 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1031 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1032 if (betar(ii) != 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1033 cnt++; |
3185 | 1034 |
10311 | 1035 ComplexColumnVector tmp (cnt); |
3185 | 1036 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1037 cnt = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1038 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1039 if (betar(ii) != 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1040 tmp(cnt++) = Complex(alphar(ii), alphai(ii))/betar(ii); |
10311 | 1041 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1042 gev = tmp; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1043 } |
3183 | 1044 } |
1045 | |
10311 | 1046 // Right, left eigenvector matrices. |
3185 | 1047 if (nargout >= 5) |
3183 | 1048 { |
10311 | 1049 // Which side to compute? |
1050 char side = (nargout == 5 ? 'R' : 'B'); | |
1051 // Compute all of them and backtransform | |
1052 char howmny = 'B'; | |
1053 // Dummy pointer; select is not used. | |
1054 octave_idx_type *select = 0; | |
3185 | 1055 |
1056 if (complex_case) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1057 { |
10311 | 1058 CVL = CQ; |
1059 CVR = CZ; | |
1060 ComplexRowVector cwork2 (2 * nn); | |
1061 RowVector rwork2 (8 * nn); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1062 octave_idx_type m; |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1063 |
10311 | 1064 F77_XFCN (ztgevc, ZTGEVC, |
1065 (F77_CONST_CHAR_ARG2 (&side, 1), | |
1066 F77_CONST_CHAR_ARG2 (&howmny, 1), | |
1067 select, nn, caa.fortran_vec (), nn, cbb.fortran_vec (), | |
1068 nn, CVL.fortran_vec (), nn, CVR.fortran_vec (), nn, nn, | |
1069 m, cwork2.fortran_vec (), rwork2.fortran_vec (), info | |
1070 F77_CHAR_ARG_LEN (1) | |
1071 F77_CHAR_ARG_LEN (1))); | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1072 } |
3185 | 1073 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1074 { |
3185 | 1075 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1076 std::cout << "qz: computing generalized eigenvectors" << std::endl; |
3185 | 1077 #endif |
1078 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1079 VL = QQ; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1080 VR = ZZ; |
10311 | 1081 octave_idx_type m; |
1082 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1083 F77_XFCN (dtgevc, DTGEVC, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1084 (F77_CONST_CHAR_ARG2 (&side, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1085 F77_CONST_CHAR_ARG2 (&howmny, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1086 select, nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1087 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
|
1088 m, work.fortran_vec (), info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1089 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1090 F77_CHAR_ARG_LEN (1))); |
3185 | 1091 |
10311 | 1092 // Now construct the complex form of VV, WW. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1093 int jj = 0; |
3185 | 1094 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1095 while (jj < nn) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1096 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1097 OCTAVE_QUIT; |
4153 | 1098 |
10311 | 1099 // See if real or complex eigenvalue. |
1100 | |
1101 // Column increment; assume complex eigenvalue. | |
1102 int cinc = 2; | |
3185 | 1103 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1104 if (jj == (nn-1)) |
10311 | 1105 // Single column. |
1106 cinc = 1; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1107 else if (aa(jj+1,jj) == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1108 cinc = 1; |
3185 | 1109 |
10311 | 1110 // Now copy the eigenvector (s) to CVR, CVL. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1111 if (cinc == 1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1112 { |
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 CVR(ii,jj) = VR(ii,jj); |
3185 | 1115 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1116 if (side == 'B') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1117 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1118 CVL(ii,jj) = VL(ii,jj); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1119 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1120 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1121 { |
10311 | 1122 // Double column; complex vector. |
3185 | 1123 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1124 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1125 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1126 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
|
1127 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
|
1128 } |
3183 | 1129 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1130 if (side == 'B') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1131 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1132 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1133 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
|
1134 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
|
1135 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1136 } |
3185 | 1137 |
10311 | 1138 // Advance to next eigenvectors (if any). |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1139 jj += cinc; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1140 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1141 } |
3183 | 1142 } |
3185 | 1143 |
1144 switch (nargout) | |
1145 { | |
1146 case 7: | |
1147 retval(6) = gev; | |
1148 | |
10311 | 1149 case 6: |
1150 // Return eigenvectors. | |
3185 | 1151 retval(5) = CVL; |
1152 | |
10311 | 1153 case 5: |
1154 // Return eigenvectors. | |
3185 | 1155 retval(4) = CVR; |
1156 | |
1157 case 4: | |
1158 if (nargin == 3) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1159 { |
3185 | 1160 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1161 std::cout << "qz: sort: retval(3) = gev = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1162 octave_print_internal (std::cout, gev); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1163 std::cout << std::endl; |
3185 | 1164 #endif |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1165 retval(3) = gev; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1166 } |
3185 | 1167 else |
10311 | 1168 { |
1169 if (complex_case) | |
1170 retval(3) = CZ; | |
1171 else | |
1172 retval(3) = ZZ; | |
1173 } | |
3185 | 1174 |
1175 case 3: | |
1176 if (nargin == 3) | |
10311 | 1177 retval(2) = CZ; |
3185 | 1178 else |
10311 | 1179 { |
1180 if (complex_case) | |
1181 retval(2) = CQ.hermitian (); | |
1182 else | |
1183 retval(2) = QQ.transpose (); | |
1184 } | |
1185 | |
3185 | 1186 case 2: |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1187 { |
10311 | 1188 if (complex_case) |
1189 { | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1190 #ifdef DEBUG |
10311 | 1191 std::cout << "qz: retval (1) = cbb = " << std::endl; |
1192 octave_print_internal (std::cout, cbb, 0); | |
1193 std::cout << std::endl << "qz: retval(0) = caa = " <<std::endl; | |
1194 octave_print_internal (std::cout, caa, 0); | |
1195 std::cout << std::endl; | |
1196 #endif | |
1197 retval(1) = cbb; | |
1198 retval(0) = caa; | |
1199 } | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1200 else |
10311 | 1201 { |
3185 | 1202 #ifdef DEBUG |
10311 | 1203 std::cout << "qz: retval (1) = bb = " << std::endl; |
1204 octave_print_internal (std::cout, bb, 0); | |
1205 std::cout << std::endl << "qz: retval(0) = aa = " <<std::endl; | |
1206 octave_print_internal (std::cout, aa, 0); | |
1207 std::cout << std::endl; | |
3185 | 1208 #endif |
10311 | 1209 retval(1) = bb; |
1210 retval(0) = aa; | |
1211 } | |
1212 } | |
1213 break; | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1214 |
3185 | 1215 |
1216 case 1: | |
1217 case 0: | |
1218 #ifdef DEBUG | |
3538 | 1219 std::cout << "qz: retval(0) = gev = " << gev << std::endl; |
3185 | 1220 #endif |
1221 retval(0) = gev; | |
1222 break; | |
1223 | |
1224 default: | |
1225 error ("qz: too many return arguments."); | |
1226 break; | |
3183 | 1227 } |
1228 | |
3185 | 1229 #ifdef DEBUG |
3538 | 1230 std::cout << "qz: exiting (at long last)" << std::endl; |
3185 | 1231 #endif |
3183 | 1232 |
1233 return retval; | |
1234 } |