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