3183
|
1 /* |
|
2 |
|
3 Copyright (C) 1998 A. S. Hodel |
|
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 |
|
9 Free Software Foundation; either version 2, or (at your option) any |
|
10 later version. |
|
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 |
|
18 along with Octave; see the file COPYING. If not, write to the Free |
|
19 Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
|
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 |
|
31 #include "config.h" |
|
32 |
|
33 #include <cfloat> |
4051
|
34 #include <cmath> |
|
35 |
3523
|
36 #include <iostream> |
4051
|
37 #include <iomanip> |
3183
|
38 |
|
39 #include "CmplxQRP.h" |
|
40 #include "dbleQR.h" |
4153
|
41 #include "f77-fcn.h" |
|
42 #include "quit.h" |
|
43 |
3183
|
44 #include "defun-dld.h" |
|
45 #include "error.h" |
|
46 #include "gripes.h" |
|
47 #include "oct-obj.h" |
|
48 #include "oct-map.h" |
|
49 #include "ov.h" |
|
50 #include "pager.h" |
3185
|
51 #if defined (DEBUG) || defined (DEBUG_SORT) |
3183
|
52 #include "pr-output.h" |
|
53 #endif |
|
54 #include "symtab.h" |
|
55 #include "utils.h" |
|
56 #include "variables.h" |
|
57 |
3185
|
58 typedef int (*sort_function) (const int& LSIZE, const double& ALPHA, |
|
59 const double& BETA, const double& S, |
|
60 const double& P); |
3183
|
61 |
|
62 extern "C" |
|
63 { |
4552
|
64 F77_RET_T |
|
65 F77_FUNC (dggbal, DGGBAL) (F77_CONST_CHAR_ARG_DECL, |
|
66 const int& N, double* A, const int& LDA, |
|
67 double* B, const int& LDB, int& ILO, |
|
68 int& IHI, double* LSCALE, double* RSCALE, |
|
69 double* WORK, int& INFO |
|
70 F77_CHAR_ARG_LEN_DECL); |
3183
|
71 |
4552
|
72 F77_RET_T |
|
73 F77_FUNC (dggbak, DGGBAK) (F77_CONST_CHAR_ARG_DECL, |
|
74 F77_CONST_CHAR_ARG_DECL, |
|
75 const int& N, const int& ILO, |
4566
|
76 const int& IHI, const double* LSCALE, |
|
77 const double* RSCALE, int& M, double* V, |
4552
|
78 const int& LDV, int& INFO |
|
79 F77_CHAR_ARG_LEN_DECL |
|
80 F77_CHAR_ARG_LEN_DECL); |
3183
|
81 |
4552
|
82 F77_RET_T |
|
83 F77_FUNC (dgghrd, DGGHRD) (F77_CONST_CHAR_ARG_DECL, |
|
84 F77_CONST_CHAR_ARG_DECL, |
|
85 const int& N, const int& ILO, |
|
86 const int& IHI, double* A, |
|
87 const int& LDA, double* B, |
|
88 const int& LDB, double* Q, |
|
89 const int& LDQ, double* Z, |
|
90 const int& LDZ, int& INFO |
|
91 F77_CHAR_ARG_LEN_DECL |
|
92 F77_CHAR_ARG_LEN_DECL); |
3183
|
93 |
4552
|
94 F77_RET_T |
|
95 F77_FUNC (dhgeqz, DHGEQZ) (F77_CONST_CHAR_ARG_DECL, |
|
96 F77_CONST_CHAR_ARG_DECL, |
|
97 F77_CONST_CHAR_ARG_DECL, |
|
98 const int& N, const int& ILO, const int& IHI, |
|
99 double* A, const int& LDA, double* B, |
|
100 const int& LDB, double* ALPHAR, |
|
101 double* ALPHAI, double* BETA, double* Q, |
|
102 const int& LDQ, double* Z, |
|
103 const int& LDZ, double* WORK, |
|
104 const int& LWORK, int& INFO |
|
105 F77_CHAR_ARG_LEN_DECL |
|
106 F77_CHAR_ARG_LEN_DECL |
|
107 F77_CHAR_ARG_LEN_DECL); |
3183
|
108 |
4552
|
109 F77_RET_T |
|
110 F77_FUNC (dlag2, DLAG2) (double* A, const int& LDA, double* B, |
|
111 const int& LDB, const double& SAFMIN, |
|
112 double& SCALE1, double& SCALE2, |
|
113 double& WR1, double& WR2, double& WI); |
3183
|
114 |
|
115 // Van Dooren's code (netlib.org: toms/590) for reordering |
|
116 // GEP. Only processes Z, not Q. |
4552
|
117 F77_RET_T |
|
118 F77_FUNC (dsubsp, DSUBSP) (const int& NMAX, const int& N, double* A, |
|
119 double* B, double* Z, sort_function, |
|
120 const double& EPS, int& NDIM, int& FAIL, |
|
121 int* IND); |
3183
|
122 |
|
123 // documentation for DTGEVC incorrectly states that VR, VL are |
|
124 // complex*16; they are declared in DTGEVC as double precision |
|
125 // (probably a cut and paste problem fro ZTGEVC) |
4552
|
126 F77_RET_T |
|
127 F77_FUNC (dtgevc, DTGEVC) (F77_CONST_CHAR_ARG_DECL, |
|
128 F77_CONST_CHAR_ARG_DECL, |
|
129 int* SELECT, const int& N, double* A, |
|
130 const int& LDA, double* B, |
|
131 const int& LDB, double* VL, |
|
132 const int& LDVL, double* VR, |
|
133 const int& LDVR, const int& MM, |
|
134 int& M, double* WORK, int& INFO |
|
135 F77_CHAR_ARG_LEN_DECL |
|
136 F77_CHAR_ARG_LEN_DECL); |
3183
|
137 |
4552
|
138 F77_RET_T |
|
139 F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG_DECL, |
|
140 double& retval |
|
141 F77_CHAR_ARG_LEN_DECL); |
3185
|
142 |
4552
|
143 F77_RET_T |
|
144 F77_FUNC (xdlange, XDLANGE) (F77_CONST_CHAR_ARG_DECL, |
|
145 const int&, const int&, const double*, |
|
146 const int&, double*, double& |
|
147 F77_CHAR_ARG_LEN_DECL); |
3183
|
148 } |
|
149 |
|
150 // fcrhp, fin, fout, folhp: |
|
151 // routines for ordering of generalized eigenvalues |
|
152 // return 1 if test is passed, 0 otherwise |
|
153 // fin: |lambda| < 1 |
|
154 // fout: |lambda| >= 1 |
|
155 // fcrhp: real(lambda) >= 0 |
|
156 // folhp: real(lambda) < 0 |
|
157 |
3185
|
158 static int |
|
159 fcrhp (const int& lsize, const double& alpha, |
|
160 const double& beta, const double& s, const double&) |
3183
|
161 { |
3185
|
162 if (lsize == 1) |
3183
|
163 return (alpha*beta >= 0 ? 1 : -1); |
3185
|
164 else |
3183
|
165 return (s >= 0 ? 1 : -1); |
|
166 } |
3185
|
167 |
|
168 static int |
|
169 fin (const int& lsize, const double& alpha, |
|
170 const double& beta, const double&, const double& p) |
3183
|
171 { |
|
172 int retval; |
|
173 |
3185
|
174 if (lsize == 1) |
|
175 retval = (fabs (alpha) < fabs (beta) ? 1 : -1); |
|
176 else |
|
177 retval = (fabs (p) < 1 ? 1 : -1); |
|
178 |
|
179 #ifdef DEBUG |
3538
|
180 std::cout << "qz: fin: retval=" << retval << std::endl; |
3185
|
181 #endif |
|
182 |
3183
|
183 return retval; |
|
184 } |
3185
|
185 |
|
186 static int |
|
187 folhp (const int& lsize, const double& alpha, |
|
188 const double& beta, const double& s, const double&) |
3183
|
189 { |
3185
|
190 if (lsize == 1) |
3183
|
191 return (alpha*beta < 0 ? 1 : -1); |
3185
|
192 else |
3183
|
193 return (s < 0 ? 1 : -1); |
|
194 } |
3185
|
195 |
|
196 static int |
|
197 fout (const int& lsize, const double& alpha, |
|
198 const double& beta, const double&, const double& p) |
3183
|
199 { |
3185
|
200 if (lsize == 1) |
|
201 return (fabs (alpha) >= fabs (beta) ? 1 : -1); |
|
202 else |
|
203 return (fabs (p) >= 1 ? 1 : -1); |
3183
|
204 } |
|
205 |
|
206 DEFUN_DLD (qz, args, nargout, |
3372
|
207 "-*- texinfo -*-\n\ |
|
208 @deftypefn {Loadable Function} {@var{lambda} =} qz (@var{a}, @var{b})\n\ |
|
209 Generalized eigenvalue problem @math{A x = s B x},\n\ |
|
210 @var{QZ} decomposition. Three ways to call:\n\ |
|
211 @enumerate\n\ |
|
212 @item @code{lambda = qz(A,B)}\n\ |
3185
|
213 \n\ |
3372
|
214 Computes the generalized eigenvalues @var{lambda} of @math{(A - sB)}.\n\ |
|
215 \n\ |
3500
|
216 @item @code{[AA, BB, Q, Z, V, W, lambda] = qz (A, B)}\n\ |
3185
|
217 \n\ |
3372
|
218 Computes qz decomposition, generalized eigenvectors, and \n\ |
|
219 generalized eigenvalues of @math{(A - sB)}\n\ |
|
220 @example\n\ |
|
221 @group\n\ |
|
222 A V = B V diag(lambda)\n\ |
|
223 W' A = diag(lambda) W' B\n\ |
|
224 AA = Q'*A*Z, BB = Q'*B*Z with Q, Z orthogonal (unitary)= I\n\ |
|
225 @end group\n\ |
|
226 @end example\n\ |
3185
|
227 \n\ |
3372
|
228 @item @code{[AA,BB,Z@{,lambda@}] = qz(A,B,opt)}\n\ |
3185
|
229 \n\ |
3372
|
230 As in form [2], but allows ordering of generalized eigenpairs\n\ |
|
231 for (e.g.) solution of discrete time algebraic Riccati equations.\n\ |
|
232 Form 3 is not available for complex matrices and does not compute\n\ |
|
233 the generalized eigenvectors V, W, nor the orthogonal matrix Q.\n\ |
|
234 @table @var\n\ |
|
235 @item opt\n\ |
|
236 for ordering eigenvalues of the GEP pencil. The leading block\n\ |
|
237 of the revised pencil contains all eigenvalues that satisfy:\n\ |
|
238 @table @code\n\ |
|
239 @item \"N\"\n\ |
|
240 = unordered (default) \n\ |
3183
|
241 \n\ |
3372
|
242 @item \"S\"\n\ |
|
243 = small: leading block has all |lambda| <=1 \n\ |
3185
|
244 \n\ |
3372
|
245 @item \"B\"\n\ |
|
246 = big: leading block has all |lambda >= 1 \n\ |
|
247 \n\ |
|
248 @item \"-\"\n\ |
|
249 = negative real part: leading block has all eigenvalues\n\ |
3183
|
250 in the open left half-plant\n\ |
3372
|
251 \n\ |
|
252 @item \"+\"\n\ |
|
253 = nonnegative real part: leading block has all eigenvalues\n\ |
3183
|
254 in the closed right half-plane\n\ |
3372
|
255 @end table\n\ |
|
256 @end table\n\ |
|
257 @end enumerate\n\ |
3183
|
258 \n\ |
3372
|
259 Note: qz performs permutation balancing, but not scaling (see balance).\n\ |
3183
|
260 Order of output arguments was selected for compatibility with MATLAB\n\ |
|
261 \n\ |
3372
|
262 See also: balance, dare, eig, schur\n\ |
|
263 @end deftypefn") |
3183
|
264 { |
|
265 octave_value_list retval; |
|
266 int nargin = args.length (); |
|
267 |
3185
|
268 #ifdef DEBUG |
3538
|
269 std::cout << "qz: nargin = " << nargin << ", nargout = " << nargout << std::endl; |
3185
|
270 #endif |
3183
|
271 |
3185
|
272 if (nargin < 2 || nargin > 3 || nargout > 7) |
|
273 { |
|
274 print_usage ("qz"); |
|
275 return retval; |
|
276 } |
|
277 else if (nargin == 3 && (nargout < 3 || nargout > 4)) |
|
278 { |
3427
|
279 error ("qz: invalid number of output arguments for form [3] call"); |
3185
|
280 return retval; |
|
281 } |
3183
|
282 |
3185
|
283 #ifdef DEBUG |
3538
|
284 std::cout << "qz: determine ordering option" << std::endl; |
3185
|
285 #endif |
3183
|
286 |
|
287 // Determine ordering option |
3523
|
288 std::string ord_job; |
3183
|
289 static double safmin; |
3185
|
290 |
|
291 if (nargin == 2) |
3183
|
292 ord_job = "N"; |
3185
|
293 else if (!args(2).is_string ()) |
|
294 { |
|
295 error ("qz: argument 3 must be a string"); |
|
296 return retval; |
|
297 } |
|
298 else |
|
299 { |
|
300 ord_job = args(2).string_value (); |
3183
|
301 |
3185
|
302 if (ord_job[0] != 'N' |
|
303 && ord_job[0] != 'S' |
|
304 && ord_job[0] != 'B' |
|
305 && ord_job[0] != '+' |
|
306 && ord_job[0] != '-') |
|
307 { |
3427
|
308 error ("qz: invalid order option"); |
3185
|
309 return retval; |
|
310 } |
|
311 |
|
312 // overflow constant required by dlag2 |
4552
|
313 F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG2 ("S", 1), |
|
314 safmin |
|
315 F77_CHAR_ARG_LEN (1)); |
3183
|
316 |
3185
|
317 #ifdef DEBUG_EIG |
3538
|
318 std::cout << "qz: initial value of safmin=" << setiosflags (std::ios::scientific) |
|
319 << safmin << std::endl; |
3185
|
320 #endif |
3183
|
321 |
3185
|
322 // some machines (e.g., DEC alpha) get safmin = 0; |
|
323 // for these, use eps instead to avoid problems in dlag2 |
|
324 if (safmin == 0) |
|
325 { |
|
326 #ifdef DEBUG_EIG |
3538
|
327 std::cout << "qz: DANGER WILL ROBINSON: safmin is 0!" << std::endl; |
3185
|
328 #endif |
3183
|
329 |
4552
|
330 F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG2 ("E", 1), |
|
331 safmin |
|
332 F77_CHAR_ARG_LEN (1)); |
3185
|
333 |
|
334 #ifdef DEBUG_EIG |
3538
|
335 std::cout << "qz: safmin set to " << setiosflags (std::ios::scientific) |
|
336 << safmin << std::endl; |
3185
|
337 #endif |
|
338 } |
3183
|
339 } |
|
340 |
3185
|
341 #ifdef DEBUG |
3538
|
342 std::cout << "qz: check argument 1" << std::endl; |
3185
|
343 #endif |
3183
|
344 |
|
345 // Argument 1: check if it's o.k. dimensioned |
3185
|
346 int nn = args(0).rows (); |
|
347 |
|
348 #ifdef DEBUG |
3531
|
349 std::cout << "argument 1 dimensions: (" << nn << "," << args(0).columns () << ")" |
3538
|
350 << std::endl; |
3185
|
351 #endif |
|
352 |
|
353 int arg_is_empty = empty_arg ("qz", nn, args(0).columns ()); |
|
354 |
3183
|
355 if (arg_is_empty < 0) |
3185
|
356 { |
|
357 gripe_empty_arg ("qz: parameter 1", 0); |
|
358 return retval; |
|
359 } |
3183
|
360 else if (arg_is_empty > 0) |
3185
|
361 { |
|
362 gripe_empty_arg ("qz: parameter 1; continuing", 0); |
|
363 return octave_value_list (2, Matrix ()); |
|
364 } |
|
365 else if (args(0).columns () != nn) |
|
366 { |
|
367 gripe_square_matrix_required ("qz"); |
|
368 return retval; |
|
369 } |
3183
|
370 |
|
371 // Argument 1: dimensions look good; get the value |
|
372 Matrix aa; |
|
373 ComplexMatrix caa; |
3185
|
374 |
|
375 if (args(0).is_complex_type ()) |
3183
|
376 caa = args(0).complex_matrix_value (); |
3185
|
377 else |
3183
|
378 aa = args(0).matrix_value (); |
3185
|
379 |
|
380 if (error_state) |
3183
|
381 return retval; |
|
382 |
3185
|
383 #ifdef DEBUG |
3538
|
384 std::cout << "qz: check argument 2" << std::endl; |
3185
|
385 #endif |
3183
|
386 |
|
387 // Extract argument 2 (bb, or cbb if complex) |
3185
|
388 if ((nn != args(1).columns ()) || (nn != args(1).rows ())) |
|
389 { |
|
390 gripe_nonconformant (); |
|
391 return retval; |
|
392 } |
|
393 |
3183
|
394 Matrix bb; |
|
395 ComplexMatrix cbb; |
3185
|
396 |
|
397 if (args(1).is_complex_type ()) |
3183
|
398 cbb = args(1).complex_matrix_value (); |
|
399 else |
|
400 bb = args(1).matrix_value (); |
3185
|
401 |
|
402 if (error_state) |
3183
|
403 return retval; |
|
404 |
|
405 // Both matrices loaded, now let's check what kind of arithmetic: |
|
406 //declared static to avoid compiler warnings about long jumps, vforks. |
3185
|
407 |
|
408 static int complex_case |
|
409 = (args(0).is_complex_type () || args(1).is_complex_type ()); |
3183
|
410 |
3185
|
411 if (nargin == 3 && complex_case) |
|
412 { |
|
413 error ("qz: cannot re-order complex qz decomposition."); |
|
414 return retval; |
|
415 } |
3183
|
416 |
|
417 // first, declare variables used in both the real and complex case |
|
418 Matrix QQ(nn,nn), ZZ(nn,nn), VR(nn,nn), VL(nn,nn); |
|
419 RowVector alphar(nn), alphai(nn), betar(nn); |
|
420 |
3185
|
421 ComplexMatrix CQ(nn,nn), CZ(nn,nn), CVR(nn,nn), CVL(nn,nn); |
3183
|
422 int ilo, ihi, info; |
3185
|
423 char compq = (nargout >= 3 ? 'V' : 'N'); |
|
424 char compz = (nargout >= 4 ? 'V' : 'N'); |
3183
|
425 |
3185
|
426 // initialize Q, Z to identity if we need either of them |
|
427 if (compq == 'V' || compz == 'V') |
|
428 for (int ii = 0; ii < nn; ii++) |
|
429 for (int jj = 0; jj < nn; jj++) |
4153
|
430 { |
|
431 OCTAVE_QUIT; |
|
432 QQ(ii,jj) = ZZ(ii,jj) = (ii == jj ? 1.0 : 0.0); |
|
433 } |
3183
|
434 |
3185
|
435 // always perform permutation balancing |
4552
|
436 const char bal_job = 'P'; |
3183
|
437 RowVector lscale(nn), rscale(nn), work(6*nn); |
|
438 |
3185
|
439 if (complex_case) |
|
440 { |
|
441 error ("Complex case not implemented yet"); |
|
442 return retval; |
|
443 } |
3183
|
444 else |
3185
|
445 { |
|
446 #ifdef DEBUG |
|
447 if (compq == 'V') |
3538
|
448 std::cout << "qz: performing balancing; QQ=" << std::endl << QQ << std::endl; |
3185
|
449 #endif |
3183
|
450 |
3185
|
451 F77_XFCN (dggbal, DGGBAL, |
4552
|
452 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
|
453 nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
|
454 nn, ilo, ihi, lscale.fortran_vec (), |
|
455 rscale.fortran_vec (), work.fortran_vec (), info |
|
456 F77_CHAR_ARG_LEN (1))); |
3185
|
457 |
|
458 if (f77_exception_encountered) |
|
459 { |
|
460 error ("unrecoverable error in qz (bal)"); |
|
461 return retval; |
|
462 } |
|
463 } |
3183
|
464 |
|
465 // Since we just want the balancing matrices, we can use dggbal |
|
466 // for both the real and complex cases; |
|
467 // left first |
3185
|
468 |
|
469 if (compq == 'V') |
|
470 { |
|
471 F77_XFCN (dggbak, DGGBAK, |
4552
|
472 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
|
473 F77_CONST_CHAR_ARG2 ("L", 1), |
4566
|
474 nn, ilo, ihi, lscale.data (), rscale.data (), |
|
475 nn, QQ.fortran_vec (), nn, info |
4552
|
476 F77_CHAR_ARG_LEN (1) |
|
477 F77_CHAR_ARG_LEN (1))); |
3183
|
478 |
3185
|
479 #ifdef DEBUG |
|
480 if (compq == 'V') |
3538
|
481 std::cout << "qz: balancing done; QQ=" << std::endl << QQ << std::endl; |
3185
|
482 #endif |
3183
|
483 |
3185
|
484 if (f77_exception_encountered) |
|
485 { |
|
486 error ("unrecoverable error in qz (bal-L)"); |
|
487 return retval; |
|
488 } |
3183
|
489 } |
|
490 |
|
491 // then right |
3185
|
492 if (compz == 'V') |
|
493 { |
4552
|
494 F77_XFCN (dggbak, DGGBAK, |
|
495 (F77_CONST_CHAR_ARG2 (&bal_job, 1), |
|
496 F77_CONST_CHAR_ARG2 ("R", 1), |
4566
|
497 nn, ilo, ihi, lscale.data (), rscale.data (), |
|
498 nn, ZZ.fortran_vec (), nn, info |
4552
|
499 F77_CHAR_ARG_LEN (1) |
|
500 F77_CHAR_ARG_LEN (1))); |
3183
|
501 |
3185
|
502 #ifdef DEBUG |
|
503 if (compz == 'V') |
3538
|
504 std::cout << "qz: balancing done; ZZ=" << std::endl << ZZ << std::endl; |
3185
|
505 #endif |
3183
|
506 |
3185
|
507 if (f77_exception_encountered) |
|
508 { |
|
509 error ("unrecoverable error in qz (bal-R)"); |
|
510 return retval; |
|
511 } |
|
512 } |
3183
|
513 |
|
514 static char qz_job; |
|
515 qz_job = (nargout < 2 ? 'E' : 'S'); |
3185
|
516 |
3183
|
517 if (complex_case) |
3185
|
518 { |
|
519 // complex case |
|
520 if (args(0).is_real_type ()) |
3587
|
521 caa = ComplexMatrix (aa); |
3185
|
522 |
|
523 if (args(1).is_real_type ()) |
3587
|
524 cbb = ComplexMatrix (bb); |
3185
|
525 |
|
526 if (compq == 'V') |
3587
|
527 CQ = ComplexMatrix (QQ); |
3185
|
528 |
|
529 if (compz == 'V') |
3587
|
530 CZ = ComplexMatrix (ZZ); |
3185
|
531 |
|
532 error ("complex case not done yet"); |
|
533 return retval; |
|
534 } |
3183
|
535 else // real matrices case |
3185
|
536 { |
|
537 #ifdef DEBUG |
3538
|
538 std::cout << "qz: peforming qr decomposition of bb" << std::endl; |
3185
|
539 #endif |
3183
|
540 |
3185
|
541 // compute the QR factorization of bb |
|
542 QR bqr (bb); |
|
543 |
|
544 #ifdef DEBUG |
3538
|
545 std::cout << "qz: qr (bb) done; now peforming qz decomposition" << std::endl; |
3185
|
546 #endif |
3183
|
547 |
3185
|
548 bb = bqr.R (); |
|
549 |
|
550 #ifdef DEBUG |
3538
|
551 std::cout << "qz: extracted bb" << std::endl; |
3185
|
552 #endif |
3183
|
553 |
3185
|
554 aa = (bqr.Q ()).transpose ()*aa; |
|
555 |
|
556 #ifdef DEBUG |
3538
|
557 std::cout << "qz: updated aa " << std::endl; |
|
558 std::cout << "bqr.Q () = " << std::endl << bqr.Q () << std::endl; |
3183
|
559 |
3185
|
560 if (compq == 'V') |
3538
|
561 std::cout << "QQ =" << QQ << std::endl; |
3185
|
562 #endif |
3183
|
563 |
3185
|
564 if (compq == 'V') |
|
565 QQ = QQ*bqr.Q (); |
3183
|
566 |
3185
|
567 #ifdef DEBUG |
3538
|
568 std::cout << "qz: precursors done..." << std::endl; |
3185
|
569 #endif |
3183
|
570 |
3185
|
571 #ifdef DEBUG |
3538
|
572 std::cout << "qz: compq = " << compq << ", compz = " << compz << std::endl; |
3185
|
573 #endif |
3183
|
574 |
3185
|
575 // reduce to generalized hessenberg form |
|
576 F77_XFCN (dgghrd, DGGHRD, |
4552
|
577 (F77_CONST_CHAR_ARG2 (&compq, 1), |
|
578 F77_CONST_CHAR_ARG2 (&compz, 1), |
|
579 nn, ilo, ihi, aa.fortran_vec (), |
|
580 nn, bb.fortran_vec (), nn, QQ.fortran_vec (), nn, |
|
581 ZZ.fortran_vec (), nn, info |
|
582 F77_CHAR_ARG_LEN (1) |
|
583 F77_CHAR_ARG_LEN (1))); |
3183
|
584 |
3185
|
585 if (f77_exception_encountered) |
|
586 { |
|
587 error ("unrecoverable error in qz (dgghrd)"); |
|
588 return retval; |
|
589 } |
3183
|
590 |
3185
|
591 // check if just computing generalized eigenvalues or if we're |
|
592 // actually computing the decomposition |
3183
|
593 |
3185
|
594 // reduce to generalized Schur form |
|
595 F77_XFCN (dhgeqz, DHGEQZ, |
4552
|
596 (F77_CONST_CHAR_ARG2 (&qz_job, 1), |
|
597 F77_CONST_CHAR_ARG2 (&compq, 1), |
|
598 F77_CONST_CHAR_ARG2 (&compz, 1), |
|
599 nn, ilo, ihi, aa.fortran_vec (), nn, bb.fortran_vec (), |
|
600 nn, alphar.fortran_vec (), alphai.fortran_vec (), |
|
601 betar.fortran_vec (), QQ.fortran_vec (), nn, |
|
602 ZZ.fortran_vec (), nn, work.fortran_vec (), nn, info |
|
603 F77_CHAR_ARG_LEN (1) |
|
604 F77_CHAR_ARG_LEN (1) |
|
605 F77_CHAR_ARG_LEN (1))); |
3183
|
606 |
3185
|
607 if (f77_exception_encountered) |
|
608 { |
|
609 error ("unrecoverable error in qz (dhgeqz)"); |
|
610 return retval; |
|
611 } |
|
612 } |
3183
|
613 |
|
614 // order the QZ decomposition? |
3185
|
615 if (ord_job[0] != 'N') |
3183
|
616 { |
3185
|
617 if (complex_case) |
|
618 { |
|
619 // probably not needed, but better be safe |
|
620 error ("qz: cannot re-order complex qz decomposition."); |
|
621 return retval; |
|
622 } |
|
623 else |
|
624 { |
|
625 #ifdef DEBUG_SORT |
3538
|
626 std::cout << "qz: ordering eigenvalues: ord_job = " << ord_job[0] << std::endl; |
3185
|
627 #endif |
3183
|
628 |
3185
|
629 // declared static to avoid vfork/long jump compiler complaints |
|
630 static sort_function sort_test; |
|
631 sort_test = NULL; |
3183
|
632 |
3185
|
633 switch (ord_job[0]) |
|
634 { |
|
635 case 'S': |
|
636 sort_test = &fin; |
|
637 break; |
3183
|
638 |
3185
|
639 case 'B': |
|
640 sort_test = &fout; |
|
641 break; |
3183
|
642 |
3185
|
643 case '+': |
|
644 sort_test = &fcrhp; |
|
645 break; |
3183
|
646 |
3185
|
647 case '-': |
|
648 sort_test = &folhp; |
|
649 break; |
|
650 |
|
651 default: |
3427
|
652 // invalid order option (should never happen, since we |
3185
|
653 // checked the options at the top). |
|
654 panic_impossible (); |
|
655 break; |
3550
|
656 } |
3183
|
657 |
3550
|
658 int ndim, fail; |
3185
|
659 double inf_norm; |
|
660 |
|
661 F77_XFCN (xdlange, XDLANGE, |
4552
|
662 (F77_CONST_CHAR_ARG2 ("I", 1), |
4566
|
663 nn, nn, aa.data (), nn, work.fortran_vec (), inf_norm |
4552
|
664 F77_CHAR_ARG_LEN (1))); |
3185
|
665 |
|
666 double eps = DBL_EPSILON*inf_norm*nn; |
|
667 |
|
668 #ifdef DEBUG_SORT |
3538
|
669 std::cout << "qz: calling dsubsp: aa=" << std::endl; |
3531
|
670 octave_print_internal (std::cout, aa, 0); |
3538
|
671 std::cout << std::endl << "bb=" << std::endl; |
3531
|
672 octave_print_internal (std::cout, bb, 0); |
3185
|
673 if (compz == 'V') |
|
674 { |
3538
|
675 std::cout << std::endl << "ZZ=" << std::endl; |
3531
|
676 octave_print_internal (std::cout, ZZ, 0); |
3185
|
677 } |
3538
|
678 std::cout << std::endl; |
|
679 std::cout << "alphar = " << std::endl; |
3531
|
680 octave_print_internal (std::cout, (Matrix) alphar, 0); |
3538
|
681 std::cout << std::endl << "alphai = " << std::endl; |
3531
|
682 octave_print_internal (std::cout, (Matrix) alphai, 0); |
3538
|
683 std::cout << std::endl << "beta = " << std::endl; |
3531
|
684 octave_print_internal (std::cout, (Matrix) betar, 0); |
3538
|
685 std::cout << std::endl; |
3185
|
686 #endif |
|
687 |
3550
|
688 Array<int> ind (nn); |
|
689 |
3185
|
690 F77_XFCN (dsubsp, DSUBSP, |
4552
|
691 (nn, nn, aa.fortran_vec (), bb.fortran_vec (), |
|
692 ZZ.fortran_vec (), sort_test, eps, ndim, fail, |
3550
|
693 ind.fortran_vec ())); |
3185
|
694 |
|
695 #ifdef DEBUG |
3538
|
696 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
3531
|
697 octave_print_internal (std::cout, aa, 0); |
3538
|
698 std::cout << std::endl << "bb=" << std::endl; |
3531
|
699 octave_print_internal (std::cout, bb, 0); |
3185
|
700 if (compz == 'V') |
|
701 { |
3538
|
702 std::cout << std::endl << "ZZ=" << std::endl; |
3531
|
703 octave_print_internal (std::cout, ZZ, 0); |
3185
|
704 } |
3538
|
705 std::cout << std::endl; |
3185
|
706 #endif |
|
707 |
|
708 // manually update alphar, alphai, betar |
|
709 static int jj; |
|
710 |
|
711 jj=0; |
|
712 while (jj < nn) |
|
713 { |
|
714 #ifdef DEBUG_EIG |
3538
|
715 std::cout << "computing gen eig #" << jj << std::endl; |
3185
|
716 #endif |
|
717 |
|
718 static int zcnt; // number of zeros in this block |
|
719 |
|
720 if (jj == (nn-1)) |
|
721 zcnt = 1; |
|
722 else if (aa(jj+1,jj) == 0) |
|
723 zcnt = 1; |
|
724 else zcnt = 2; |
|
725 |
|
726 if (zcnt == 1) // real zero |
|
727 { |
|
728 #ifdef DEBUG_EIG |
3538
|
729 std::cout << " single gen eig:" << std::endl; |
|
730 std::cout << " alphar(" << jj << ") = " << aa(jj,jj) << std::endl; |
|
731 std::cout << " betar( " << jj << ") = " << bb(jj,jj) << std::endl; |
|
732 std::cout << " alphai(" << jj << ") = 0" << std::endl; |
3185
|
733 #endif |
|
734 |
|
735 alphar(jj) = aa(jj,jj); |
|
736 alphai(jj) = 0; |
|
737 betar(jj) = bb(jj,jj); |
|
738 } |
|
739 else |
|
740 { |
|
741 // complex conjugate pair |
|
742 #ifdef DEBUG_EIG |
3538
|
743 std::cout << "qz: calling dlag2:" << std::endl; |
3531
|
744 std::cout << "safmin=" |
3538
|
745 << setiosflags (std::ios::scientific) << safmin << std::endl; |
3185
|
746 |
|
747 for (int idr = jj; idr <= jj+1; idr++) |
|
748 { |
|
749 for (int idc = jj; idc <= jj+1; idc++) |
|
750 { |
3531
|
751 std::cout << "aa(" << idr << "," << idc << ")=" |
3538
|
752 << aa(idr,idc) << std::endl; |
3531
|
753 std::cout << "bb(" << idr << "," << idc << ")=" |
3538
|
754 << bb(idr,idc) << std::endl; |
3185
|
755 } |
|
756 } |
|
757 #endif |
|
758 |
4566
|
759 // XXX FIXME XXX -- probably should be using |
|
760 // fortran_vec instead of &aa(jj,jj) here. |
|
761 |
3185
|
762 double scale1, scale2, wr1, wr2, wi; |
|
763 F77_XFCN (dlag2, DLAG2, |
|
764 (&aa(jj,jj), nn, &bb(jj,jj), nn, safmin, |
|
765 scale1, scale2, wr1, wr2, wi)); |
|
766 |
|
767 #ifdef DEBUG_EIG |
3531
|
768 std::cout << "dlag2 returns: scale1=" << scale1 |
3538
|
769 << "\tscale2=" << scale2 << std::endl |
3185
|
770 << "\twr1=" << wr1 << "\twr2=" << wr2 |
3538
|
771 << "\twi=" << wi << std::endl; |
3185
|
772 #endif |
|
773 |
|
774 // just to be safe, check if it's a real pair |
|
775 if (wi == 0) |
|
776 { |
|
777 alphar(jj) = wr1; |
|
778 alphai(jj) = 0; |
|
779 betar(jj) = scale1; |
|
780 alphar(jj+1) = wr2; |
|
781 alphai(jj+1) = 0; |
|
782 betar(jj+1) = scale2; |
|
783 } |
|
784 else |
|
785 { |
|
786 alphar(jj) = alphar(jj+1)=wr1; |
|
787 alphai(jj) = -(alphai(jj+1) = wi); |
|
788 betar(jj) = betar(jj+1) = scale1; |
|
789 } |
|
790 } |
|
791 |
|
792 // advance past this block |
|
793 jj += zcnt; |
|
794 } |
|
795 |
|
796 #ifdef DEBUG_SORT |
3538
|
797 std::cout << "qz: back from dsubsp: aa=" << std::endl; |
3531
|
798 octave_print_internal (std::cout, aa, 0); |
3538
|
799 std::cout << std::endl << "bb=" << std::endl; |
3531
|
800 octave_print_internal (std::cout, bb, 0); |
3185
|
801 |
|
802 if (compz == 'V') |
|
803 { |
3538
|
804 std::cout << std::endl << "ZZ=" << std::endl; |
3531
|
805 octave_print_internal (std::cout, ZZ, 0); |
3185
|
806 } |
3538
|
807 std::cout << std::endl << "qz: ndim=" << ndim << std::endl |
|
808 << "fail=" << fail << std::endl; |
|
809 std::cout << "alphar = " << std::endl; |
3531
|
810 octave_print_internal (std::cout, (Matrix) alphar, 0); |
3538
|
811 std::cout << std::endl << "alphai = " << std::endl; |
3531
|
812 octave_print_internal (std::cout, (Matrix) alphai, 0); |
3538
|
813 std::cout << std::endl << "beta = " << std::endl; |
3531
|
814 octave_print_internal (std::cout, (Matrix) betar, 0); |
3538
|
815 std::cout << std::endl; |
3185
|
816 #endif |
|
817 } |
3183
|
818 } |
3185
|
819 |
3183
|
820 // compute generalized eigenvalues? |
|
821 ComplexColumnVector gev; |
3185
|
822 |
|
823 if (nargout < 2 || nargout == 7 || (nargin == 3 && nargout == 4)) |
3183
|
824 { |
3185
|
825 if (complex_case) |
|
826 { |
|
827 error ("complex case not yet implemented"); |
|
828 return retval; |
|
829 } |
|
830 else |
|
831 { |
|
832 #ifdef DEBUG |
3538
|
833 std::cout << "qz: computing generalized eigenvalues" << std::endl; |
3185
|
834 #endif |
3183
|
835 |
3185
|
836 // return finite generalized eigenvalues |
|
837 int cnt = 0; |
|
838 |
|
839 for (int ii = 0; ii < nn; ii++) |
|
840 if (betar(ii) != 0) |
|
841 cnt++; |
|
842 |
|
843 ComplexColumnVector tmp(cnt); |
|
844 |
3671
|
845 cnt = 0; |
3185
|
846 for (int ii = 0; ii < nn; ii++) |
|
847 if (betar(ii) != 0) |
3671
|
848 tmp(cnt++) = Complex(alphar(ii), alphai(ii))/betar(ii); |
3185
|
849 gev = tmp; |
|
850 } |
3183
|
851 } |
|
852 |
|
853 // right, left eigenvector matrices |
3185
|
854 if (nargout >= 5) |
3183
|
855 { |
3185
|
856 char side = (nargout == 5 ? 'R' : 'B'); // which side to compute? |
|
857 char howmny = 'B'; // compute all of them and backtransform |
|
858 int *select = NULL; // dummy pointer; select is not used. |
|
859 int m; |
|
860 |
|
861 if (complex_case) |
|
862 { |
|
863 error ("complex type not yet implemented"); |
|
864 return retval; |
|
865 } |
|
866 else |
|
867 { |
|
868 #ifdef DEBUG |
3538
|
869 std::cout << "qz: computing generalized eigenvectors" << std::endl; |
3185
|
870 #endif |
|
871 |
|
872 VL = QQ; |
|
873 VR = ZZ; |
|
874 |
|
875 F77_XFCN (dtgevc, DTGEVC, |
4552
|
876 (F77_CONST_CHAR_ARG2 (&side, 1), |
|
877 F77_CONST_CHAR_ARG2 (&howmny, 1), |
|
878 select, nn, aa.fortran_vec (), nn, bb.fortran_vec (), |
|
879 nn, VL.fortran_vec (), nn, VR.fortran_vec (), nn, nn, |
|
880 m, work.fortran_vec (), info |
|
881 F77_CHAR_ARG_LEN (1) |
|
882 F77_CHAR_ARG_LEN (1))); |
3185
|
883 |
|
884 if (f77_exception_encountered) |
|
885 { |
|
886 error ("unrecoverable error in qz (dtgevc)"); |
|
887 return retval; |
|
888 } |
3183
|
889 |
3185
|
890 // now construct the complex form of VV, WW |
|
891 int jj = 0; |
|
892 |
|
893 while (jj < nn) |
|
894 { |
4153
|
895 OCTAVE_QUIT; |
|
896 |
3185
|
897 // see if real or complex eigenvalue |
|
898 int cinc = 2; // column increment; assume complex eigenvalue |
|
899 |
|
900 if (jj == (nn-1)) |
|
901 cinc = 1; // single column |
|
902 else if (aa(jj+1,jj) == 0) |
|
903 cinc = 1; |
|
904 |
|
905 // now copy the eigenvector (s) to CVR, CVL |
|
906 if (cinc == 1) |
|
907 { |
|
908 for (int ii = 0; ii < nn; ii++) |
|
909 CVR(ii,jj) = VR(ii,jj); |
|
910 |
|
911 if (side == 'B') |
|
912 for (int ii = 0; ii < nn; ii++) |
|
913 CVL(ii,jj) = VL(ii,jj); |
|
914 } |
|
915 else |
|
916 { |
|
917 // double column; complex vector |
|
918 |
|
919 for (int ii = 0; ii < nn; ii++) |
|
920 { |
|
921 CVR(ii,jj) = Complex (VR(ii,jj), VR(ii,jj+1)); |
|
922 CVR(ii,jj+1) = Complex (VR(ii,jj), -VR(ii,jj+1)); |
|
923 } |
3183
|
924 |
3185
|
925 if (side == 'B') |
|
926 for (int ii = 0; ii < nn; ii++) |
|
927 { |
|
928 CVL(ii,jj) = Complex (VL(ii,jj), VL(ii,jj+1)); |
|
929 CVL(ii,jj+1) = Complex (VL(ii,jj), -VL(ii,jj+1)); |
|
930 } |
|
931 } |
|
932 |
|
933 // advance to next eigenvectors (if any) |
|
934 jj += cinc; |
|
935 } |
|
936 } |
3183
|
937 } |
3185
|
938 |
|
939 switch (nargout) |
|
940 { |
|
941 case 7: |
|
942 retval(6) = gev; |
|
943 |
|
944 case 6: // return eigenvectors |
|
945 retval(5) = CVL; |
|
946 |
|
947 case 5: // return eigenvectors |
|
948 retval(4) = CVR; |
|
949 |
|
950 case 4: |
|
951 if (nargin == 3) |
|
952 { |
|
953 #ifdef DEBUG |
3538
|
954 std::cout << "qz: sort: retval(3) = gev = " << std::endl; |
3531
|
955 octave_print_internal (std::cout, gev); |
3538
|
956 std::cout << std::endl; |
3185
|
957 #endif |
|
958 retval(3) = gev; |
|
959 } |
|
960 else |
|
961 retval(3) = ZZ; |
|
962 |
|
963 case 3: |
|
964 if (nargin == 3) |
|
965 retval(2) = ZZ; |
|
966 else |
|
967 retval(2) = QQ; |
|
968 |
|
969 case 2: |
|
970 #ifdef DEBUG |
3538
|
971 std::cout << "qz: retval (1) = bb = " << std::endl; |
3531
|
972 octave_print_internal (std::cout, bb, 0); |
3538
|
973 std::cout << std::endl << "qz: retval(0) = aa = " <<std::endl; |
3531
|
974 octave_print_internal (std::cout, aa, 0); |
3538
|
975 std::cout << std::endl; |
3185
|
976 #endif |
|
977 retval(1) = bb; |
|
978 retval(0) = aa; |
|
979 break; |
|
980 |
|
981 case 1: |
|
982 case 0: |
|
983 #ifdef DEBUG |
3538
|
984 std::cout << "qz: retval(0) = gev = " << gev << std::endl; |
3185
|
985 #endif |
|
986 retval(0) = gev; |
|
987 break; |
|
988 |
|
989 default: |
|
990 error ("qz: too many return arguments."); |
|
991 break; |
3183
|
992 } |
|
993 |
3185
|
994 #ifdef DEBUG |
3538
|
995 std::cout << "qz: exiting (at long last)" << std::endl; |
3185
|
996 #endif |
3183
|
997 |
|
998 return retval; |
|
999 } |
|
1000 |
|
1001 /* |
|
1002 ;;; Local Variables: *** |
|
1003 ;;; mode: C++ *** |
|
1004 ;;; End: *** |
|
1005 */ |