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