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