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