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