Mercurial > octave-nkf
annotate libinterp/corefcn/qz.cc @ 19437:03067dab10ca
Use stricter input validation when looking for a string as input (bug #42651).
* data.cc (get_sort_mode_option, Fissorted): Use is_string() to check string
input.
* debug.cc (Fdbstep): use "string" rather than "character string" in error
messages.
* error.cc (Flasterr, Flastwarn): use "string" rather than "character string"
in error messages.
* file-io.cc (do_stream_open, do_fread, do_fwrite, Fpopen, Ftempname,
Fmkstemp): Use is_string() to check string input.
* graphics.cc (Fgraphics_toolkit): Use is_string() to check string input.
Rephrase error message.
* help.cc (F__list_functions): Use is_string() to check string input.
* input.cc (Fyes_or_no): Use is_string() to check string input. Rephrase
error message.
* input.cc (Fadd_input_event_hook): Rephrase error message.
* load-path.cc (Fgenpath, Faddpath): Rephrase error message.
* matrix_type.cc (Fmatrix_type): Use is_string() to check string input.
* qz.cc (Fqz): Follow Octave coding convention for space after '!'.
* regexp.cc (parse_options): Use is_string() to check string input.
Rephrase error message.
* schur.cc (Fschur): Use is_string() to check string input.
* strfns.cc (Flist_in_columns): Use is_string() to check string input.
Rephrase error message.
* symtab.cc (Fignore_function_time_stamp): Use is_string() to check string
input. Rephrase error message.
* syscalls.cc (Fexec, Fpopen2, Fcanonicalize_file_name): Use is_string() to
check string input. Rephrase error message.
* sysdep.cc (Fsetenv): Use is_string() to check string input.
* time.cc (Fstrftime, Fstrptime): Use is_string() to check string input.
* toplev.cc (Fsystem, Fatexit): Use is_string() to check string input.
* urlwrite.cc (Furlwrite, Furlread): Rephrase error message.
* utils.cc (Ffile_in_path): Use is_string() to check string input. Rephrase
error message.
* variables.cc (extract_function): Add FIXME about potentially using is_string.
* variables.cc (do_isglobal, Fmunlock, Fmislocked): Use is_string() to check
string input.
* variables.cc (set_internal_variable): Rephrase error message.
* ov-base.cc (make_idx_args): Rephrase error message.
* ov-class.cc (octave_class::all_strings, Fclass): Rephrase error message.
* ov-fcn-handle.cc (Fstr2func): Use is_string() to check string input
* ov-java.cc (FjavaObject, FjavaMethod, F__java_get__, F__java_set__):
Use is_string() to check string input.
* ov.cc (Fdecode_subscripts): Use is_string() to check string input.
Rephrase error message.
* pt-idx.cc (tree_index_expression::get_struct_index): Rephrase error message.
* io.tst: Change %!warning test to %!error test to match stricter checking.
* system.tst: Change %!warning test for setenv to %!error test to match
stricter checking.
author | Rik <rik@octave.org> |
---|---|
date | Tue, 16 Dec 2014 09:21:29 -0800 |
parents | 6113e0c6920b |
children | 4197fc428c7d |
rev | line source |
---|---|
3183 | 1 /* |
2 | |
17744
d63878346099
maint: Update copyright notices for release.
John W. Eaton <jwe@octave.org>
parents:
17281
diff
changeset
|
3 Copyright (C) 1998-2013 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 |
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 | |
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 |
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, | |
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 |
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 | |
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\ |
15039
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14876
diff
changeset
|
296 @deftypefn {Built-in Function} {@var{lambda} =} qz (@var{A}, @var{B})\n\ |
e753177cde93
maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents:
14876
diff
changeset
|
297 @deftypefnx {Built-in Function} {@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\ |
1577c6f80926
Use non-breaking spaces between certain adjectives and their nouns in docstrings.
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
299 (@math{A x = s B x}). There are three ways to call this function:\n\ |
3372 | 300 @enumerate\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
301 @item @code{@var{lambda} = qz (@var{A}, @var{B})}\n\ |
3185 | 302 \n\ |
5016 | 303 Computes the generalized eigenvalues\n\ |
304 @tex\n\ | |
305 $\\lambda$\n\ | |
306 @end tex\n\ | |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
307 @ifnottex\n\ |
5016 | 308 @var{lambda}\n\ |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
309 @end ifnottex\n\ |
5016 | 310 of @math{(A - s B)}.\n\ |
10840 | 311 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
312 @item @code{[AA, BB, Q, Z, V, W, @var{lambda}] = qz (@var{A}, @var{B})}\n\ |
3185 | 313 \n\ |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
11593
diff
changeset
|
314 Computes QZ@tie{}decomposition, generalized eigenvectors, and\n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
315 generalized eigenvalues of @math{(A - s B)}\n\ |
5016 | 316 @tex\n\ |
317 $$ AV = BV{ \\rm diag }(\\lambda) $$\n\ | |
318 $$ W^T A = { \\rm diag }(\\lambda)W^T B $$\n\ | |
319 $$ AA = Q^T AZ, BB = Q^T BZ $$\n\ | |
320 @end tex\n\ | |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
321 @ifnottex\n\ |
10840 | 322 \n\ |
3372 | 323 @example\n\ |
324 @group\n\ | |
5481 | 325 \n\ |
14360
97883071e8e4
doc: Correct off-by-1 spacings in all .cc docstrings
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
326 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
|
327 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
|
328 AA = Q * A * Z, BB = Q * B * Z\n\ |
5481 | 329 \n\ |
3372 | 330 @end group\n\ |
331 @end example\n\ | |
10840 | 332 \n\ |
8517
81d6ab3ac93c
Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents:
8286
diff
changeset
|
333 @end ifnottex\n\ |
5016 | 334 with @var{Q} and @var{Z} orthogonal (unitary)= @var{I}\n\ |
3185 | 335 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
336 @item @code{[AA,BB,Z@{, @var{lambda}@}] = qz (@var{A}, @var{B}, @var{opt})}\n\ |
3185 | 337 \n\ |
3372 | 338 As in form [2], but allows ordering of generalized eigenpairs\n\ |
5481 | 339 for (e.g.) solution of discrete time algebraic Riccati equations.\n\ |
340 Form 3 is not available for complex matrices, and does not compute\n\ | |
10840 | 341 the generalized eigenvectors @var{V}, @var{W}, nor the orthogonal matrix\n\ |
342 @var{Q}.\n\ | |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10350
diff
changeset
|
343 \n\ |
3372 | 344 @table @var\n\ |
345 @item opt\n\ | |
16826
a4969508008e
doc: Periodic spellcheck of the documentation.
Rik <rik@octave.org>
parents:
16772
diff
changeset
|
346 for ordering eigenvalues of the @nospell{GEP} pencil. The leading block\n\ |
5481 | 347 of 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
|
348 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
349 @table @asis\n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
350 @item @qcode{\"N\"}\n\ |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
11593
diff
changeset
|
351 = unordered (default)\n\ |
3183 | 352 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
353 @item @qcode{\"S\"}\n\ |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
11593
diff
changeset
|
354 = small: leading block has all |lambda| @leq{} 1\n\ |
3185 | 355 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
356 @item @qcode{\"B\"}\n\ |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
11593
diff
changeset
|
357 = big: leading block has all |lambda| @geq{} 1\n\ |
3372 | 358 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
359 @item @qcode{\"-\"}\n\ |
5481 | 360 = negative real part: leading block has all eigenvalues\n\ |
361 in the open left half-plane\n\ | |
3372 | 362 \n\ |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
363 @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
|
364 = non-negative real part: leading block has all eigenvalues\n\ |
5481 | 365 in the closed right half-plane\n\ |
9064
7c02ec148a3c
Check grammar on all .cc files
Rik <rdrider0-list@yahoo.com>
parents:
8927
diff
changeset
|
366 @end table\n\ |
3372 | 367 @end table\n\ |
368 @end enumerate\n\ | |
3183 | 369 \n\ |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
370 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
|
371 (@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
|
372 compatibility with @sc{matlab}.\n\ |
16920
53eaa83e4181
doc: Add seealso links between various factorization forms.
Rik <rik@octave.org>
parents:
16826
diff
changeset
|
373 @seealso{eig, balance, lu, chol, hess, qr, qzhess, schur, svd}\n\ |
3372 | 374 @end deftypefn") |
3183 | 375 { |
376 octave_value_list retval; | |
377 int nargin = args.length (); | |
378 | |
3185 | 379 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
380 std::cout << "qz: nargin = " << nargin |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
381 << ", nargout = " << nargout << std::endl; |
3185 | 382 #endif |
3183 | 383 |
3185 | 384 if (nargin < 2 || nargin > 3 || nargout > 7) |
385 { | |
5823 | 386 print_usage (); |
3185 | 387 return retval; |
388 } | |
389 else if (nargin == 3 && (nargout < 3 || nargout > 4)) | |
390 { | |
3427 | 391 error ("qz: invalid number of output arguments for form [3] call"); |
3185 | 392 return retval; |
393 } | |
3183 | 394 |
3185 | 395 #ifdef DEBUG |
3538 | 396 std::cout << "qz: determine ordering option" << std::endl; |
3185 | 397 #endif |
3183 | 398 |
10311 | 399 // 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
|
400 volatile char ord_job = 0; |
3183 | 401 static double safmin; |
3185 | 402 |
403 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
|
404 ord_job = 'N'; |
19437
03067dab10ca
Use stricter input validation when looking for a string as input (bug #42651).
Rik <rik@octave.org>
parents:
18712
diff
changeset
|
405 else if (! args(2).is_string ()) |
3185 | 406 { |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
407 error ("qz: OPT must be a string"); |
3185 | 408 return retval; |
409 } | |
410 else | |
411 { | |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
412 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
|
413 |
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
414 if (! tmp.empty ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
415 ord_job = tmp[0]; |
3183 | 416 |
8927
d75f4ee0538d
qz.cc (Fqz): avoid maybe clobbred by vfork warning from GCC
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
417 if (! (ord_job == 'N' || ord_job == 'n' |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
418 || ord_job == 'S' || ord_job == 's' |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
419 || ord_job == 'B' || ord_job == 'b' |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
420 || ord_job == '+' || ord_job == '-')) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
421 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
422 error ("qz: invalid order option"); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
423 return retval; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
424 } |
3185 | 425 |
426 // overflow constant required by dlag2 | |
4552 | 427 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
|
428 safmin |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
429 F77_CHAR_ARG_LEN (1)); |
3183 | 430 |
3185 | 431 #ifdef DEBUG_EIG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
432 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
|
433 << setiosflags (std::ios::scientific) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
434 << safmin << std::endl; |
3185 | 435 #endif |
3183 | 436 |
10311 | 437 // Some machines (e.g., DEC alpha) get safmin = 0; |
438 // for these, use eps instead to avoid problems in dlag2. | |
3185 | 439 if (safmin == 0) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
440 { |
3185 | 441 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
442 std::cout << "qz: DANGER WILL ROBINSON: safmin is 0!" << std::endl; |
3185 | 443 #endif |
3183 | 444 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
445 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
|
446 safmin |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
447 F77_CHAR_ARG_LEN (1)); |
3185 | 448 |
449 #ifdef DEBUG_EIG | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
450 std::cout << "qz: safmin set to " |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
451 << setiosflags (std::ios::scientific) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
452 << safmin << std::endl; |
3185 | 453 #endif |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
454 } |
3183 | 455 } |
456 | |
3185 | 457 #ifdef DEBUG |
3538 | 458 std::cout << "qz: check argument 1" << std::endl; |
3185 | 459 #endif |
3183 | 460 |
10311 | 461 // Argument 1: check if it's o.k. dimensioned. |
5275 | 462 octave_idx_type nn = args(0).rows (); |
3185 | 463 |
464 #ifdef DEBUG | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
465 std::cout << "argument 1 dimensions: (" |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
466 << nn << "," << args(0).columns () << ")" |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
467 << std::endl; |
3185 | 468 #endif |
469 | |
470 int arg_is_empty = empty_arg ("qz", nn, args(0).columns ()); | |
471 | |
3183 | 472 if (arg_is_empty < 0) |
3185 | 473 { |
474 gripe_empty_arg ("qz: parameter 1", 0); | |
475 return retval; | |
476 } | |
3183 | 477 else if (arg_is_empty > 0) |
3185 | 478 { |
479 gripe_empty_arg ("qz: parameter 1; continuing", 0); | |
480 return octave_value_list (2, Matrix ()); | |
481 } | |
482 else if (args(0).columns () != nn) | |
483 { | |
484 gripe_square_matrix_required ("qz"); | |
485 return retval; | |
486 } | |
3183 | 487 |
10311 | 488 // Argument 1: dimensions look good; get the value. |
3183 | 489 Matrix aa; |
490 ComplexMatrix caa; | |
3185 | 491 |
492 if (args(0).is_complex_type ()) | |
3183 | 493 caa = args(0).complex_matrix_value (); |
3185 | 494 else |
3183 | 495 aa = args(0).matrix_value (); |
3185 | 496 |
497 if (error_state) | |
3183 | 498 return retval; |
499 | |
3185 | 500 #ifdef DEBUG |
3538 | 501 std::cout << "qz: check argument 2" << std::endl; |
3185 | 502 #endif |
3183 | 503 |
10311 | 504 // Extract argument 2 (bb, or cbb if complex). |
3185 | 505 if ((nn != args(1).columns ()) || (nn != args(1).rows ())) |
506 { | |
507 gripe_nonconformant (); | |
508 return retval; | |
509 } | |
510 | |
3183 | 511 Matrix bb; |
512 ComplexMatrix cbb; | |
3185 | 513 |
514 if (args(1).is_complex_type ()) | |
3183 | 515 cbb = args(1).complex_matrix_value (); |
516 else | |
517 bb = args(1).matrix_value (); | |
3185 | 518 |
519 if (error_state) | |
3183 | 520 return retval; |
521 | |
522 // Both matrices loaded, now let's check what kind of arithmetic: | |
10311 | 523 // declared volatile to avoid compiler warnings about long jumps, |
524 // vforks. | |
3185 | 525 |
10311 | 526 volatile int complex_case |
3185 | 527 = (args(0).is_complex_type () || args(1).is_complex_type ()); |
10311 | 528 |
3185 | 529 if (nargin == 3 && complex_case) |
530 { | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
531 error ("qz: cannot re-order complex qz decomposition"); |
3185 | 532 return retval; |
533 } | |
3183 | 534 |
10311 | 535 // First, declare variables used in both the real and complex case. |
3183 | 536 Matrix QQ(nn,nn), ZZ(nn,nn), VR(nn,nn), VL(nn,nn); |
537 RowVector alphar(nn), alphai(nn), betar(nn); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
538 ComplexRowVector xalpha(nn), xbeta(nn); |
3185 | 539 ComplexMatrix CQ(nn,nn), CZ(nn,nn), CVR(nn,nn), CVL(nn,nn); |
5275 | 540 octave_idx_type ilo, ihi, info; |
3185 | 541 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
|
542 char compz = ((nargout >= 4 || nargin == 3)? 'V' : 'N'); |
3183 | 543 |
10311 | 544 // Initialize Q, Z to identity if we need either of them. |
3185 | 545 if (compq == 'V' || compz == 'V') |
5275 | 546 for (octave_idx_type ii = 0; ii < nn; ii++) |
547 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
|
548 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
549 OCTAVE_QUIT; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
550 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
|
551 } |
3183 | 552 |
10311 | 553 // Always perform permutation balancing. |
4552 | 554 const char bal_job = 'P'; |
10311 | 555 RowVector lscale (nn), rscale (nn), work (6 * nn), rwork (nn); |
3183 | 556 |
3185 | 557 if (complex_case) |
558 { | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
559 #ifdef DEBUG |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
560 if (compq == 'V') |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
561 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
|
562 << CQ << std::endl; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
563 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
564 if (args(0).is_real_type ()) |
10311 | 565 caa = ComplexMatrix (aa); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
566 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
567 if (args(1).is_real_type ()) |
10311 | 568 cbb = ComplexMatrix (bb); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
569 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
570 if (compq == 'V') |
10311 | 571 CQ = ComplexMatrix (QQ); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
572 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
573 if (compz == 'V') |
10311 | 574 CZ = ComplexMatrix (ZZ); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
575 |
10311 | 576 F77_XFCN (zggbal, ZGGBAL, |
577 (F77_CONST_CHAR_ARG2 (&bal_job, 1), | |
578 nn, caa.fortran_vec (), nn, cbb.fortran_vec (), | |
579 nn, ilo, ihi, lscale.fortran_vec (), | |
580 rscale.fortran_vec (), work.fortran_vec (), info | |
581 F77_CHAR_ARG_LEN (1))); | |
3185 | 582 } |
3183 | 583 else |
3185 | 584 { |
585 #ifdef DEBUG | |
586 if (compq == 'V') | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
587 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
|
588 << QQ << std::endl; |
3185 | 589 #endif |
3183 | 590 |
3185 | 591 F77_XFCN (dggbal, DGGBAL, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
592 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
593 nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
594 nn, ilo, ihi, lscale.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
595 rscale.fortran_vec (), work.fortran_vec (), info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
596 F77_CHAR_ARG_LEN (1))); |
3185 | 597 } |
3183 | 598 |
599 // Since we just want the balancing matrices, we can use dggbal | |
10311 | 600 // for both the real and complex cases; left first |
3185 | 601 |
10311 | 602 #if 0 |
603 if (compq == 'V') | |
3185 | 604 { |
605 F77_XFCN (dggbak, DGGBAK, | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
606 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
607 F77_CONST_CHAR_ARG2 ("L", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
608 nn, ilo, ihi, lscale.data (), rscale.data (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
609 nn, QQ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
610 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
611 F77_CHAR_ARG_LEN (1))); |
3183 | 612 |
3185 | 613 #ifdef DEBUG |
614 if (compq == 'V') | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
615 std::cout << "qz: balancing done; QQ=" << std::endl << QQ << std::endl; |
3185 | 616 #endif |
3183 | 617 } |
618 | |
619 // then right | |
3185 | 620 if (compz == 'V') |
621 { | |
4552 | 622 F77_XFCN (dggbak, DGGBAK, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
623 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
624 F77_CONST_CHAR_ARG2 ("R", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
625 nn, ilo, ihi, lscale.data (), rscale.data (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
626 nn, ZZ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
627 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
628 F77_CHAR_ARG_LEN (1))); |
3183 | 629 |
3185 | 630 #ifdef DEBUG |
631 if (compz == 'V') | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
632 std::cout << "qz: balancing done; ZZ=" << std::endl << ZZ << std::endl; |
3185 | 633 #endif |
10311 | 634 } |
635 #endif | |
3183 | 636 |
637 static char qz_job; | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
638 qz_job = (nargout < 2 ? 'E' : 'S'); |
3185 | 639 |
3183 | 640 if (complex_case) |
3185 | 641 { |
10311 | 642 // Complex case. |
643 | |
644 // The QR decomposition of cbb. | |
645 ComplexQR cbqr (cbb); | |
646 // The R matrix of QR decomposition for cbb. | |
647 cbb = cbqr.R (); | |
648 // (Q*)caa for following work. | |
649 caa = (cbqr.Q ().hermitian ()) * caa; | |
650 CQ = CQ * cbqr.Q (); | |
651 | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
652 F77_XFCN (zgghrd, ZGGHRD, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
653 (F77_CONST_CHAR_ARG2 (&compq, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
654 F77_CONST_CHAR_ARG2 (&compz, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
655 nn, ilo, ihi, caa.fortran_vec (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
656 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
|
657 CZ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
658 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
659 F77_CHAR_ARG_LEN (1))); |
10311 | 660 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
661 ComplexRowVector cwork (1 * nn); |
10311 | 662 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
663 F77_XFCN (zhgeqz, ZHGEQZ, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
664 (F77_CONST_CHAR_ARG2 (&qz_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
665 F77_CONST_CHAR_ARG2 (&compq, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
666 F77_CONST_CHAR_ARG2 (&compz, 1), |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
667 nn, ilo, ihi, |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
668 caa.fortran_vec (), nn, |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
669 cbb.fortran_vec (),nn, |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
670 xalpha.fortran_vec (), xbeta.fortran_vec (), |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
671 CQ.fortran_vec (), nn, |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
672 CZ.fortran_vec (), nn, |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
673 cwork.fortran_vec (), nn, rwork.fortran_vec (), info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
674 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
675 F77_CHAR_ARG_LEN (1) |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
676 F77_CHAR_ARG_LEN (1))); |
3185 | 677 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
678 if (compq == 'V') |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
679 { |
10311 | 680 // Left eigenvector. |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
681 F77_XFCN (zggbak, ZGGBAK, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
682 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
683 F77_CONST_CHAR_ARG2 ("L", 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
684 nn, ilo, ihi, lscale.data (), rscale.data (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
685 nn, CQ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
686 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
687 F77_CHAR_ARG_LEN (1))); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
688 } |
3185 | 689 |
10311 | 690 // Right eigenvector. |
3185 | 691 if (compz == 'V') |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
692 { |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
693 F77_XFCN (zggbak, ZGGBAK, |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
694 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
695 F77_CONST_CHAR_ARG2 ("R", 1), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
696 nn, ilo, ihi, lscale.data (), rscale.data (), |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
697 nn, CZ.fortran_vec (), nn, info |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
698 F77_CHAR_ARG_LEN (1) |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
699 F77_CHAR_ARG_LEN (1))); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
700 } |
3185 | 701 |
702 } | |
10311 | 703 else |
3185 | 704 { |
705 #ifdef DEBUG | |
3538 | 706 std::cout << "qz: peforming qr decomposition of bb" << std::endl; |
3185 | 707 #endif |
3183 | 708 |
10311 | 709 // Compute the QR factorization of bb. |
3185 | 710 QR bqr (bb); |
711 | |
712 #ifdef DEBUG | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
713 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
|
714 << std::endl; |
3185 | 715 #endif |
3183 | 716 |
3185 | 717 bb = bqr.R (); |
718 | |
719 #ifdef DEBUG | |
3538 | 720 std::cout << "qz: extracted bb" << std::endl; |
3185 | 721 #endif |
3183 | 722 |
10311 | 723 aa = (bqr.Q ()).transpose () * aa; |
3185 | 724 |
725 #ifdef DEBUG | |
3538 | 726 std::cout << "qz: updated aa " << std::endl; |
727 std::cout << "bqr.Q () = " << std::endl << bqr.Q () << std::endl; | |
3183 | 728 |
3185 | 729 if (compq == 'V') |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
730 std::cout << "QQ =" << QQ << std::endl; |
3185 | 731 #endif |
3183 | 732 |
3185 | 733 if (compq == 'V') |
10311 | 734 QQ = QQ * bqr.Q (); |
3183 | 735 |
3185 | 736 #ifdef DEBUG |
3538 | 737 std::cout << "qz: precursors done..." << std::endl; |
3185 | 738 #endif |
3183 | 739 |
3185 | 740 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
741 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
|
742 << std::endl; |
3185 | 743 #endif |
3183 | 744 |
10311 | 745 // Reduce to generalized hessenberg form. |
3185 | 746 F77_XFCN (dgghrd, DGGHRD, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
747 (F77_CONST_CHAR_ARG2 (&compq, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
748 F77_CONST_CHAR_ARG2 (&compz, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
749 nn, ilo, ihi, aa.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
750 nn, bb.fortran_vec (), nn, QQ.fortran_vec (), nn, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
751 ZZ.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
752 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
753 F77_CHAR_ARG_LEN (1))); |
3183 | 754 |
10311 | 755 // Check if just computing generalized eigenvalues or if we're |
756 // actually computing the decomposition. | |
3183 | 757 |
10311 | 758 // Reduce to generalized Schur form. |
3185 | 759 F77_XFCN (dhgeqz, DHGEQZ, |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
760 (F77_CONST_CHAR_ARG2 (&qz_job, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
761 F77_CONST_CHAR_ARG2 (&compq, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
762 F77_CONST_CHAR_ARG2 (&compz, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
763 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
|
764 nn, alphar.fortran_vec (), alphai.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
765 betar.fortran_vec (), QQ.fortran_vec (), nn, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
766 ZZ.fortran_vec (), nn, work.fortran_vec (), nn, info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
767 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
768 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
769 F77_CHAR_ARG_LEN (1))); |
10311 | 770 |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
771 if (compq == 'V') |
10311 | 772 { |
773 F77_XFCN (dggbak, DGGBAK, | |
774 (F77_CONST_CHAR_ARG2 (&bal_job, 1), | |
775 F77_CONST_CHAR_ARG2 ("L", 1), | |
776 nn, ilo, ihi, lscale.data (), rscale.data (), | |
777 nn, QQ.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 |
10311 | 782 if (compq == 'V') |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
783 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
|
784 << QQ << std::endl; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
785 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
786 } |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
787 |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
788 // then right |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
789 if (compz == 'V') |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
790 { |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
791 F77_XFCN (dggbak, DGGBAK, |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
792 (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
|
793 F77_CONST_CHAR_ARG2 ("R", 1), |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
794 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
|
795 nn, ZZ.fortran_vec (), nn, info |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
796 F77_CHAR_ARG_LEN (1) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
797 F77_CHAR_ARG_LEN (1))); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
798 |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
799 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
800 if (compz == 'V') |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
801 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
|
802 << ZZ << std::endl; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
803 #endif |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
804 } |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
805 |
3185 | 806 } |
3183 | 807 |
10311 | 808 // 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
|
809 if (! (ord_job == 'N' || ord_job == 'n')) |
3183 | 810 { |
3185 | 811 if (complex_case) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
812 { |
10311 | 813 // Probably not needed, but better be safe. |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
814 error ("qz: cannot re-order complex qz decomposition"); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
815 return retval; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
816 } |
3185 | 817 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
818 { |
3185 | 819 #ifdef DEBUG_SORT |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
820 std::cout << "qz: ordering eigenvalues: ord_job = " |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
821 << ord_job << std::endl; |
3185 | 822 #endif |
3183 | 823 |
10311 | 824 // 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
|
825 static sort_function sort_test; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
826 sort_test = 0; |
3183 | 827 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
828 switch (ord_job) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
829 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
830 case 'S': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
831 case 's': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
832 sort_test = &fin; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
833 break; |
3183 | 834 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
835 case 'B': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
836 case 'b': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
837 sort_test = &fout; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
838 break; |
3183 | 839 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
840 case '+': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
841 sort_test = &fcrhp; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
842 break; |
3183 | 843 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
844 case '-': |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
845 sort_test = &folhp; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
846 break; |
3185 | 847 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
848 default: |
10311 | 849 // Invalid order option (should never happen, since we |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
850 // checked the options at the top). |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
851 panic_impossible (); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
852 break; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
853 } |
3183 | 854 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
855 octave_idx_type ndim, fail; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
856 double inf_norm; |
3185 | 857 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
858 F77_XFCN (xdlange, XDLANGE, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
859 (F77_CONST_CHAR_ARG2 ("I", 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
860 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
|
861 F77_CHAR_ARG_LEN (1))); |
3185 | 862 |
15220
61822c866ba1
use std::numeric_limits<T>::epsilon in C++ code
John W. Eaton <jwe@octave.org>
parents:
15195
diff
changeset
|
863 double eps = std::numeric_limits<double>::epsilon () * inf_norm * nn; |
3185 | 864 |
865 #ifdef DEBUG_SORT | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
866 std::cout << "qz: calling dsubsp: aa=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
867 octave_print_internal (std::cout, aa, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
868 std::cout << std::endl << "bb=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
869 octave_print_internal (std::cout, bb, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
870 if (compz == 'V') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
871 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
872 std::cout << std::endl << "ZZ=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
873 octave_print_internal (std::cout, ZZ, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
874 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
875 std::cout << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
876 std::cout << "alphar = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
877 octave_print_internal (std::cout, (Matrix) alphar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
878 std::cout << std::endl << "alphai = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
879 octave_print_internal (std::cout, (Matrix) alphai, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
880 std::cout << std::endl << "beta = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
881 octave_print_internal (std::cout, (Matrix) betar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
882 std::cout << std::endl; |
3185 | 883 #endif |
884 | |
11570
57632dea2446
attempt better backward compatibility for Array constructors
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
885 Array<octave_idx_type> ind (dim_vector (nn, 1)); |
3550 | 886 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
887 F77_XFCN (dsubsp, DSUBSP, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
888 (nn, nn, aa.fortran_vec (), bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
889 ZZ.fortran_vec (), sort_test, eps, ndim, fail, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
890 ind.fortran_vec ())); |
3185 | 891 |
892 #ifdef DEBUG | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
893 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
894 octave_print_internal (std::cout, aa, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
895 std::cout << std::endl << "bb=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
896 octave_print_internal (std::cout, bb, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
897 if (compz == 'V') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
898 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
899 std::cout << std::endl << "ZZ=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
900 octave_print_internal (std::cout, ZZ, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
901 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
902 std::cout << std::endl; |
3185 | 903 #endif |
904 | |
10311 | 905 // Manually update alphar, alphai, betar. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
906 static int jj; |
3185 | 907 |
10311 | 908 jj = 0; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
909 while (jj < nn) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
910 { |
3185 | 911 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
912 std::cout << "computing gen eig #" << jj << std::endl; |
3185 | 913 #endif |
914 | |
10311 | 915 // Number of zeros in this block. |
916 static int zcnt; | |
3185 | 917 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
918 if (jj == (nn-1)) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
919 zcnt = 1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
920 else if (aa(jj+1,jj) == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
921 zcnt = 1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
922 else zcnt = 2; |
3185 | 923 |
10311 | 924 if (zcnt == 1) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
925 { |
10311 | 926 // Real zero. |
3185 | 927 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
928 std::cout << " single gen eig:" << std::endl; |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
929 std::cout << " alphar(" << jj << ") = " << aa(jj,jj) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
930 << std::endl; |
18712
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
931 std::cout << " betar(" << jj << ") = " << bb(jj,jj) |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
932 << std::endl; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
933 std::cout << " alphai(" << jj << ") = 0" << std::endl; |
3185 | 934 #endif |
935 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
936 alphar(jj) = aa(jj,jj); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
937 alphai(jj) = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
938 betar(jj) = bb(jj,jj); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
939 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
940 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
941 { |
10311 | 942 // Complex conjugate pair. |
3185 | 943 #ifdef DEBUG_EIG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
944 std::cout << "qz: calling dlag2:" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
945 std::cout << "safmin=" |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
946 << setiosflags (std::ios::scientific) |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
947 << safmin << std::endl; |
3185 | 948 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
949 for (int idr = jj; idr <= jj+1; idr++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
950 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
951 for (int idc = jj; idc <= jj+1; idc++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
952 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
953 std::cout << "aa(" << idr << "," << idc << ")=" |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
954 << aa(idr,idc) << std::endl; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
955 std::cout << "bb(" << idr << "," << idc << ")=" |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
956 << bb(idr,idc) << std::endl; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
957 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
958 } |
3185 | 959 #endif |
960 | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
961 // FIXME: probably should be using |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
962 // fortran_vec instead of &aa(jj,jj) here. |
4566 | 963 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
964 double scale1, scale2, wr1, wr2, wi; |
10311 | 965 const double *aa_ptr = aa.data () + jj * nn + jj; |
966 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
|
967 F77_XFCN (dlag2, DLAG2, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
968 (aa_ptr, nn, bb_ptr, nn, safmin, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
969 scale1, scale2, wr1, wr2, wi)); |
3185 | 970 |
971 #ifdef DEBUG_EIG | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
972 std::cout << "dlag2 returns: scale1=" << scale1 |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
973 << "\tscale2=" << scale2 << std::endl |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
974 << "\twr1=" << wr1 << "\twr2=" << wr2 |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
975 << "\twi=" << wi << std::endl; |
3185 | 976 #endif |
977 | |
10311 | 978 // 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
|
979 if (wi == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
980 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
981 alphar(jj) = wr1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
982 alphai(jj) = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
983 betar(jj) = scale1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
984 alphar(jj+1) = wr2; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
985 alphai(jj+1) = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
986 betar(jj+1) = scale2; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
987 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
988 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
989 { |
10311 | 990 alphar(jj) = alphar(jj+1) = wr1; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
991 alphai(jj) = -(alphai(jj+1) = wi); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
992 betar(jj) = betar(jj+1) = scale1; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
993 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
994 } |
3185 | 995 |
10311 | 996 // Advance past this block. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
997 jj += zcnt; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
998 } |
3185 | 999 |
1000 #ifdef DEBUG_SORT | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1001 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1002 octave_print_internal (std::cout, aa, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1003 std::cout << std::endl << "bb=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1004 octave_print_internal (std::cout, bb, 0); |
3185 | 1005 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1006 if (compz == 'V') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1007 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1008 std::cout << std::endl << "ZZ=" << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1009 octave_print_internal (std::cout, ZZ, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1010 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1011 std::cout << std::endl << "qz: ndim=" << ndim << std::endl |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1012 << "fail=" << fail << std::endl; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1013 std::cout << "alphar = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1014 octave_print_internal (std::cout, (Matrix) alphar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1015 std::cout << std::endl << "alphai = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1016 octave_print_internal (std::cout, (Matrix) alphai, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1017 std::cout << std::endl << "beta = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1018 octave_print_internal (std::cout, (Matrix) betar, 0); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1019 std::cout << std::endl; |
3185 | 1020 #endif |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1021 } |
3183 | 1022 } |
3185 | 1023 |
10311 | 1024 // Compute generalized eigenvalues? |
3183 | 1025 ComplexColumnVector gev; |
3185 | 1026 |
1027 if (nargout < 2 || nargout == 7 || (nargin == 3 && nargout == 4)) | |
3183 | 1028 { |
3185 | 1029 if (complex_case) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1030 { |
10311 | 1031 int cnt = 0; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1032 |
10311 | 1033 for (int ii = 0; ii < nn; ii++) |
1034 cnt++; | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1035 |
10311 | 1036 ComplexColumnVector tmp (cnt); |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1037 |
10311 | 1038 cnt = 0; |
1039 for (int ii = 0; ii < nn; ii++) | |
1040 tmp(cnt++) = xalpha(ii) / xbeta(ii); | |
1041 | |
1042 gev = tmp; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1043 } |
3185 | 1044 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1045 { |
3185 | 1046 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1047 std::cout << "qz: computing generalized eigenvalues" << std::endl; |
3185 | 1048 #endif |
3183 | 1049 |
10311 | 1050 // Return finite generalized eigenvalues. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1051 int cnt = 0; |
3185 | 1052 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1053 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1054 if (betar(ii) != 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1055 cnt++; |
3185 | 1056 |
10311 | 1057 ComplexColumnVector tmp (cnt); |
3185 | 1058 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1059 cnt = 0; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1060 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1061 if (betar(ii) != 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1062 tmp(cnt++) = Complex(alphar(ii), alphai(ii))/betar(ii); |
10311 | 1063 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1064 gev = tmp; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1065 } |
3183 | 1066 } |
1067 | |
10311 | 1068 // Right, left eigenvector matrices. |
3185 | 1069 if (nargout >= 5) |
3183 | 1070 { |
10311 | 1071 // Which side to compute? |
1072 char side = (nargout == 5 ? 'R' : 'B'); | |
1073 // Compute all of them and backtransform | |
1074 char howmny = 'B'; | |
1075 // Dummy pointer; select is not used. | |
1076 octave_idx_type *select = 0; | |
3185 | 1077 |
1078 if (complex_case) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1079 { |
10311 | 1080 CVL = CQ; |
1081 CVR = CZ; | |
1082 ComplexRowVector cwork2 (2 * nn); | |
1083 RowVector rwork2 (8 * nn); | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1084 octave_idx_type m; |
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1085 |
10311 | 1086 F77_XFCN (ztgevc, ZTGEVC, |
1087 (F77_CONST_CHAR_ARG2 (&side, 1), | |
1088 F77_CONST_CHAR_ARG2 (&howmny, 1), | |
1089 select, nn, caa.fortran_vec (), nn, cbb.fortran_vec (), | |
1090 nn, CVL.fortran_vec (), nn, CVR.fortran_vec (), nn, nn, | |
1091 m, cwork2.fortran_vec (), rwork2.fortran_vec (), info | |
1092 F77_CHAR_ARG_LEN (1) | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
1093 F77_CHAR_ARG_LEN (1))); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1094 } |
3185 | 1095 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1096 { |
3185 | 1097 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1098 std::cout << "qz: computing generalized eigenvectors" << std::endl; |
3185 | 1099 #endif |
1100 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1101 VL = QQ; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1102 VR = ZZ; |
10311 | 1103 octave_idx_type m; |
1104 | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1105 F77_XFCN (dtgevc, DTGEVC, |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1106 (F77_CONST_CHAR_ARG2 (&side, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1107 F77_CONST_CHAR_ARG2 (&howmny, 1), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1108 select, nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1109 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
|
1110 m, work.fortran_vec (), info |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1111 F77_CHAR_ARG_LEN (1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1112 F77_CHAR_ARG_LEN (1))); |
3185 | 1113 |
10311 | 1114 // Now construct the complex form of VV, WW. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1115 int jj = 0; |
3185 | 1116 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1117 while (jj < nn) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1118 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1119 OCTAVE_QUIT; |
4153 | 1120 |
10311 | 1121 // See if real or complex eigenvalue. |
1122 | |
1123 // Column increment; assume complex eigenvalue. | |
1124 int cinc = 2; | |
3185 | 1125 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1126 if (jj == (nn-1)) |
10311 | 1127 // Single column. |
1128 cinc = 1; | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1129 else if (aa(jj+1,jj) == 0) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1130 cinc = 1; |
3185 | 1131 |
10311 | 1132 // Now copy the eigenvector (s) to CVR, CVL. |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1133 if (cinc == 1) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1134 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1135 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1136 CVR(ii,jj) = VR(ii,jj); |
3185 | 1137 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1138 if (side == 'B') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1139 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1140 CVL(ii,jj) = VL(ii,jj); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1141 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1142 else |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1143 { |
10311 | 1144 // Double column; complex vector. |
3185 | 1145 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1146 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1147 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1148 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
|
1149 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
|
1150 } |
3183 | 1151 |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1152 if (side == 'B') |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1153 for (int ii = 0; ii < nn; ii++) |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1154 { |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1155 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
|
1156 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
|
1157 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1158 } |
3185 | 1159 |
10311 | 1160 // Advance to next eigenvectors (if any). |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1161 jj += cinc; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1162 } |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1163 } |
3183 | 1164 } |
3185 | 1165 |
1166 switch (nargout) | |
1167 { | |
1168 case 7: | |
1169 retval(6) = gev; | |
1170 | |
10311 | 1171 case 6: |
1172 // Return eigenvectors. | |
3185 | 1173 retval(5) = CVL; |
1174 | |
10311 | 1175 case 5: |
1176 // Return eigenvectors. | |
3185 | 1177 retval(4) = CVR; |
1178 | |
1179 case 4: | |
1180 if (nargin == 3) | |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1181 { |
3185 | 1182 #ifdef DEBUG |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1183 std::cout << "qz: sort: retval(3) = gev = " << std::endl; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1184 octave_print_internal (std::cout, gev); |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1185 std::cout << std::endl; |
3185 | 1186 #endif |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1187 retval(3) = gev; |
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
1188 } |
3185 | 1189 else |
10311 | 1190 { |
1191 if (complex_case) | |
1192 retval(3) = CZ; | |
1193 else | |
1194 retval(3) = ZZ; | |
1195 } | |
3185 | 1196 |
1197 case 3: | |
1198 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
|
1199 { |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1200 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
|
1201 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
|
1202 else |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1203 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
|
1204 } |
3185 | 1205 else |
10311 | 1206 { |
1207 if (complex_case) | |
1208 retval(2) = CQ.hermitian (); | |
1209 else | |
1210 retval(2) = QQ.transpose (); | |
1211 } | |
1212 | |
3185 | 1213 case 2: |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1214 { |
10311 | 1215 if (complex_case) |
1216 { | |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1217 #ifdef DEBUG |
14844
5bc9b9cb4362
maint: Use Octave coding conventions for cuddled parenthesis in retval assignments.
Rik <octave@nomad.inbox5.com>
parents:
14501
diff
changeset
|
1218 std::cout << "qz: retval(1) = cbb = " << std::endl; |
10311 | 1219 octave_print_internal (std::cout, cbb, 0); |
1220 std::cout << std::endl << "qz: retval(0) = caa = " <<std::endl; | |
1221 octave_print_internal (std::cout, caa, 0); | |
1222 std::cout << std::endl; | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
1223 #endif |
10311 | 1224 retval(1) = cbb; |
1225 retval(0) = caa; | |
1226 } | |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1227 else |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1228 { |
3185 | 1229 #ifdef DEBUG |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1230 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
|
1231 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
|
1232 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
|
1233 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
|
1234 std::cout << std::endl; |
3185 | 1235 #endif |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1236 retval(1) = bb; |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1237 retval(0) = aa; |
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1238 } |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11570
diff
changeset
|
1239 } |
10311 | 1240 break; |
10307
4e4270ab70d6
support complex case in qz
Jyh-miin Lin <jyhmiinlin@gmail.com>
parents:
10155
diff
changeset
|
1241 |
3185 | 1242 |
1243 case 1: | |
1244 case 0: | |
1245 #ifdef DEBUG | |
3538 | 1246 std::cout << "qz: retval(0) = gev = " << gev << std::endl; |
3185 | 1247 #endif |
1248 retval(0) = gev; | |
1249 break; | |
1250 | |
1251 default: | |
11553
01f703952eff
Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
1252 error ("qz: too many return arguments"); |
3185 | 1253 break; |
3183 | 1254 } |
1255 | |
3185 | 1256 #ifdef DEBUG |
3538 | 1257 std::cout << "qz: exiting (at long last)" << std::endl; |
3185 | 1258 #endif |
3183 | 1259 |
1260 return retval; | |
1261 } | |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1262 |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1263 /* |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1264 %!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
|
1265 %! 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
|
1266 %! 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
|
1267 %! 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
|
1268 %!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
|
1269 %!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
|
1270 |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1271 ## 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
|
1272 %!test |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1273 %! a = [ 10 1 2; |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1274 %! 1 2 -1; |
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1275 %! 1 1 2]; |
14501
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1276 %! 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
|
1277 %! [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
|
1278 %! sz = length (lambda); |
60e5cf354d80
Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents:
14360
diff
changeset
|
1279 %! observed = (b * v * diag ([lambda;0])) (:, 1:sz); |
18712
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
1280 %! 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
|
1281 %! observed = (diag ([lambda;0]) * w' * b) (1:sz, :); |
18712
6113e0c6920b
maint: Clean up extra spaces before/after parentheses.
Rik <rik@octave.org>
parents:
17787
diff
changeset
|
1282 %! 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
|
1283 %! 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
|
1284 %! 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
|
1285 |
14876
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1286 %!test |
54a386f2ac4e
correcly compute Z for 3-output call to qz (bug #36728)
John W. Eaton <jwe@octave.org>
parents:
14844
diff
changeset
|
1287 %! 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
|
1288 %! 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
|
1289 %! [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
|
1290 %! [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
|
1291 %! assert (Z1, Z2); |
13074
0b789a03bde1
codesprint: Add 3 tests for qz.cc
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents:
12642
diff
changeset
|
1292 */ |