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