GNU Octave  4.4.1
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
qz.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 1998-2018 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
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11 
12 Octave is distributed in the hope that it will be useful, but
13 WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License 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, see
19 <https://www.gnu.org/licenses/>.
20 
21 */
22 
23 // Generalized eigenvalue balancing via LAPACK
24 
25 // Originally written by A. S. Hodel <scotte@eng.auburn.edu>, but is
26 // substantially different with the change to use LAPACK.
27 
28 #undef DEBUG
29 #undef DEBUG_SORT
30 #undef DEBUG_EIG
31 
32 #if defined (HAVE_CONFIG_H)
33 # include "config.h"
34 #endif
35 
36 #include <cctype>
37 #include <cmath>
38 
39 #if defined (DEBUG) || defined (DEBUG_SORT) || defined (DEBUG_EIG)
40 # include <iostream>
41 # if defined (DEBUG_EIG)
42 # include <iomanip>
43 # endif
44 #endif
45 
46 #include "f77-fcn.h"
47 #include "lo-lapack-proto.h"
48 #include "qr.h"
49 #include "quit.h"
50 
51 #include "defun.h"
52 #include "error.h"
53 #include "errwarn.h"
54 #include "ovl.h"
55 #if defined (DEBUG) || defined (DEBUG_SORT)
56 # include "pr-output.h"
57 #endif
58 
59 // FIXME: Matlab does not produce lambda as the first output argument.
60 // Compatibility problem?
61 
62 DEFUN (qz, args, nargout,
63  doc: /* -*- texinfo -*-
64 @deftypefn {} {@var{lambda} =} qz (@var{A}, @var{B})
65 @deftypefnx {} {[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}, @var{lambda}] =} qz (@var{A}, @var{B})
66 @deftypefnx {} {[@var{AA}, @var{BB}, @var{Z}] =} qz (@var{A}, @var{B}, @var{opt})
67 @deftypefnx {} {[@var{AA}, @var{BB}, @var{Z}, @var{lambda}] =} qz (@var{A}, @var{B}, @var{opt})
68 Compute the QZ@tie{}decomposition of a generalized eigenvalue problem.
69 
70 The generalized eigenvalue problem is defined as
71 
72 @tex
73 $$A x = \lambda B x$$
74 @end tex
75 @ifnottex
76 
77 @math{A x = @var{lambda} B x}
78 
79 @end ifnottex
80 
81 There are three calling forms of the function:
82 
83 @enumerate
84 @item @code{@var{lambda} = qz (@var{A}, @var{B})}
85 
86 Compute the generalized eigenvalues
87 @tex
88 $\lambda.$
89 @end tex
90 @ifnottex
91 @var{lambda}.
92 @end ifnottex
93 
94 @item @code{[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}, @var{lambda}] = qz (@var{A}, @var{B})}
95 
96 Compute QZ@tie{}decomposition, generalized eigenvectors, and generalized
97 eigenvalues.
98 @tex
99 $$ AV = BV{ \rm diag }(\lambda) $$
100 $$ W^T A = { \rm diag }(\lambda)W^T B $$
101 $$ AA = Q^T AZ, BB = Q^T BZ $$
102 @end tex
103 @ifnottex
104 
105 @example
106 @group
107 
108 @var{A} * @var{V} = @var{B} * @var{V} * diag (@var{lambda})
109 @var{W}' * @var{A} = diag (@var{lambda}) * @var{W}' * @var{B}
110 @var{AA} = @var{Q} * @var{A} * @var{Z}, @var{BB} = @var{Q} * @var{B} * @var{Z}
111 
112 @end group
113 @end example
114 
115 @end ifnottex
116 with @var{Q} and @var{Z} orthogonal (unitary for complex case).
117 
118 @item @code{[@var{AA}, @var{BB}, @var{Z} @{, @var{lambda}@}] = qz (@var{A}, @var{B}, @var{opt})}
119 
120 As in form 2 above, but allows ordering of generalized eigenpairs for, e.g.,
121 solution of discrete time algebraic @nospell{Riccati} equations. Form 3 is not
122 available for complex matrices, and does not compute the generalized
123 eigenvectors @var{V}, @var{W}, nor the orthogonal matrix @var{Q}.
124 
125 @table @var
126 @item opt
127 for ordering eigenvalues of the @nospell{GEP} pencil. The leading block of
128 the revised pencil contains all eigenvalues that satisfy:
129 
130 @table @asis
131 @item @qcode{"N"}
132 unordered (default)
133 
134 @item @qcode{"S"}
135 small: leading block has all
136 @tex
137 $|\lambda| < 1$
138 @end tex
139 @ifnottex
140 |@var{lambda}| < 1
141 @end ifnottex
142 
143 @item @qcode{"B"}
144 big: leading block has all
145 @tex
146 $|\lambda| \geq 1$
147 @end tex
148 @ifnottex
149 |@var{lambda}| @geq{} 1
150 @end ifnottex
151 
152 @item @qcode{"-"}
153 negative real part: leading block has all eigenvalues in the open left
154 half-plane
155 
156 @item @qcode{"+"}
157 non-negative real part: leading block has all eigenvalues in the closed right
158 half-plane
159 @end table
160 @end table
161 @end enumerate
162 
163 Note: @code{qz} performs permutation balancing, but not scaling
164 (@pxref{XREFbalance,,balance}), which may be lead to less accurate results than
165 @code{eig}. The order of output arguments was selected for compatibility with
166 @sc{matlab}.
167 @seealso{eig, balance, lu, chol, hess, qr, qzhess, schur, svd}
168 @end deftypefn */)
169 {
170  int nargin = args.length ();
171 
172 #if defined (DEBUG)
173  std::cout << "qz: nargin = " << nargin
174  << ", nargout = " << nargout << std::endl;
175 #endif
176 
177  if (nargin < 2 || nargin > 3 || nargout > 7)
178  print_usage ();
179 
180  if (nargin == 3 && (nargout < 3 || nargout > 4))
181  error ("qz: invalid number of output arguments for form 3 call");
182 
183 #if defined (DEBUG)
184  std::cout << "qz: determine ordering option" << std::endl;
185 #endif
186 
187  // Determine ordering option.
188 
189  char ord_job = 'N';
190  double safmin = 0.0;
191 
192  if (nargin == 3)
193  {
194  std::string opt = args(2).xstring_value ("qz: OPT must be a string");
195 
196  if (opt.empty ())
197  error ("qz: OPT must be a non-empty string");
198 
199  ord_job = std::toupper (opt[0]);
200 
201  std::string valid_opts = "NSB+-";
202 
203  if (valid_opts.find_first_of (ord_job) == std::string::npos)
204  error ("qz: invalid order option '%c'", ord_job);
205 
206  // overflow constant required by dlag2
207  F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG2 ("S", 1),
208  safmin
209  F77_CHAR_ARG_LEN (1));
210 
211 #if defined (DEBUG_EIG)
212  std::cout << "qz: initial value of safmin="
213  << setiosflags (std::ios::scientific)
214  << safmin << std::endl;
215 #endif
216 
217  // Some machines (e.g., DEC alpha) get safmin = 0;
218  // for these, use eps instead to avoid problems in dlag2.
219  if (safmin == 0)
220  {
221 #if defined (DEBUG_EIG)
222  std::cout << "qz: DANGER WILL ROBINSON: safmin is 0!" << std::endl;
223 #endif
224 
225  F77_FUNC (xdlamch, XDLAMCH) (F77_CONST_CHAR_ARG2 ("E", 1),
226  safmin
227  F77_CHAR_ARG_LEN (1));
228 
229 #if defined (DEBUG_EIG)
230  std::cout << "qz: safmin set to "
231  << setiosflags (std::ios::scientific)
232  << safmin << std::endl;
233 #endif
234  }
235  }
236 
237 #if defined (DEBUG)
238  std::cout << "qz: check matrix A" << std::endl;
239 #endif
240 
241  // Matrix A: check dimensions.
242  F77_INT nn = octave::to_f77_int (args(0).rows ());
243  F77_INT nc = octave::to_f77_int (args(0).columns ());
244 
245 #if defined (DEBUG)
246  std::cout << "Matrix A dimensions: (" << nn << ',' << nc << ')' << std::endl;
247 #endif
248 
249  if (args(0).isempty ())
250  {
251  warn_empty_arg ("qz: A");
252  return octave_value_list (2, Matrix ());
253  }
254  else if (nc != nn)
255  err_square_matrix_required ("qz", "A");
256 
257  // Matrix A: get value.
258  Matrix aa;
259  ComplexMatrix caa;
260 
261  if (args(0).iscomplex ())
262  caa = args(0).complex_matrix_value ();
263  else
264  aa = args(0).matrix_value ();
265 
266 #if defined (DEBUG)
267  std::cout << "qz: check matrix B" << std::endl;
268 #endif
269 
270  // Extract argument 2 (bb, or cbb if complex).
271  F77_INT b_nr = octave::to_f77_int (args(1).rows ());
272  F77_INT b_nc = octave::to_f77_int (args(1).columns ());
273 
274  if (nn != b_nc || nn != b_nr)
276 
277  Matrix bb;
278  ComplexMatrix cbb;
279 
280  if (args(1).iscomplex ())
281  cbb = args(1).complex_matrix_value ();
282  else
283  bb = args(1).matrix_value ();
284 
285  // Both matrices loaded, now check whether to calculate complex or real val.
286 
287  bool complex_case = (args(0).iscomplex () || args(1).iscomplex ());
288 
289  if (nargin == 3 && complex_case)
290  error ("qz: cannot re-order complex qz decomposition");
291 
292  // First, declare variables used in both the real and complex cases.
293  // FIXME: There are a lot of excess variables declared.
294  // Probably a better way to handle this.
295  Matrix QQ (nn,nn), ZZ (nn,nn), VR (nn,nn), VL (nn,nn);
296  RowVector alphar (nn), alphai (nn), betar (nn);
297  ComplexRowVector xalpha (nn), xbeta (nn);
298  ComplexMatrix CQ (nn,nn), CZ (nn,nn), CVR (nn,nn), CVL (nn,nn);
299  F77_INT ilo, ihi, info;
300  char comp_q = (nargout >= 3 ? 'V' : 'N');
301  char comp_z = ((nargout >= 4 || nargin == 3)? 'V' : 'N');
302 
303  // Initialize Q, Z to identity matrix if either is needed
304  if (comp_q == 'V' || comp_z == 'V')
305  {
306  double *QQptr = QQ.fortran_vec ();
307  double *ZZptr = ZZ.fortran_vec ();
308  std::fill_n (QQptr, QQ.numel (), 0.0);
309  std::fill_n (ZZptr, ZZ.numel (), 0.0);
310  for (F77_INT i = 0; i < nn; i++)
311  {
312  QQ(i,i) = 1.0;
313  ZZ(i,i) = 1.0;
314  }
315  }
316 
317  // Always perform permutation balancing.
318  const char bal_job = 'P';
319  RowVector lscale (nn), rscale (nn), work (6 * nn), rwork (nn);
320 
321  if (complex_case)
322  {
323 #if defined (DEBUG)
324  if (comp_q == 'V')
325  std::cout << "qz: performing balancing; CQ=" << std::endl
326  << CQ << std::endl;
327 #endif
328  if (args(0).isreal ())
329  caa = ComplexMatrix (aa);
330 
331  if (args(1).isreal ())
332  cbb = ComplexMatrix (bb);
333 
334  if (comp_q == 'V')
335  CQ = ComplexMatrix (QQ);
336 
337  if (comp_z == 'V')
338  CZ = ComplexMatrix (ZZ);
339 
340  F77_XFCN (zggbal, ZGGBAL,
341  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
342  nn, F77_DBLE_CMPLX_ARG (caa.fortran_vec ()), nn,
344  nn, ilo, ihi, lscale.fortran_vec (),
345  rscale.fortran_vec (), work.fortran_vec (), info
346  F77_CHAR_ARG_LEN (1)));
347  }
348  else
349  {
350 #if defined (DEBUG)
351  if (comp_q == 'V')
352  std::cout << "qz: performing balancing; QQ=" << std::endl
353  << QQ << std::endl;
354 #endif
355 
356  F77_XFCN (dggbal, DGGBAL,
357  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
358  nn, aa.fortran_vec (), nn, bb.fortran_vec (),
359  nn, ilo, ihi, lscale.fortran_vec (),
360  rscale.fortran_vec (), work.fortran_vec (), info
361  F77_CHAR_ARG_LEN (1)));
362  }
363 
364  // Only permutation balance above is done. Skip scaling balance.
365 
366 #if 0
367  // Since we just want the balancing matrices, we can use dggbal
368  // for both the real and complex cases; left first
369 
370  if (comp_q == 'V')
371  {
372  F77_XFCN (dggbak, DGGBAK,
373  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
374  F77_CONST_CHAR_ARG2 ("L", 1),
375  nn, ilo, ihi, lscale.data (), rscale.data (),
376  nn, QQ.fortran_vec (), nn, info
377  F77_CHAR_ARG_LEN (1)
378  F77_CHAR_ARG_LEN (1)));
379 
380 #if defined (DEBUG)
381  if (comp_q == 'V')
382  std::cout << "qz: balancing done; QQ=" << std::endl << QQ << std::endl;
383 #endif
384  }
385 
386  // then right
387  if (comp_z == 'V')
388  {
389  F77_XFCN (dggbak, DGGBAK,
390  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
391  F77_CONST_CHAR_ARG2 ("R", 1),
392  nn, ilo, ihi, lscale.data (), rscale.data (),
393  nn, ZZ.fortran_vec (), nn, info
394  F77_CHAR_ARG_LEN (1)
395  F77_CHAR_ARG_LEN (1)));
396 
397 #if defined (DEBUG)
398  if (comp_z == 'V')
399  std::cout << "qz: balancing done; ZZ=" << std::endl << ZZ << std::endl;
400 #endif
401  }
402 #endif
403 
404  char qz_job = (nargout < 2 ? 'E' : 'S');
405 
406  if (complex_case)
407  {
408  // Complex case.
409 
410  // The QR decomposition of cbb.
412  // The R matrix of QR decomposition for cbb.
413  cbb = cbqr.R ();
414  // (Q*)caa for following work.
415  caa = (cbqr.Q ().hermitian ()) * caa;
416  CQ = CQ * cbqr.Q ();
417 
418  F77_XFCN (zgghrd, ZGGHRD,
419  (F77_CONST_CHAR_ARG2 (&comp_q, 1),
420  F77_CONST_CHAR_ARG2 (&comp_z, 1),
421  nn, ilo, ihi, F77_DBLE_CMPLX_ARG (caa.fortran_vec ()),
422  nn, F77_DBLE_CMPLX_ARG (cbb.fortran_vec ()), nn,
423  F77_DBLE_CMPLX_ARG (CQ.fortran_vec ()), nn,
424  F77_DBLE_CMPLX_ARG (CZ.fortran_vec ()), nn, info
425  F77_CHAR_ARG_LEN (1)
426  F77_CHAR_ARG_LEN (1)));
427 
428  ComplexRowVector cwork (nn);
429 
430  F77_XFCN (zhgeqz, ZHGEQZ,
431  (F77_CONST_CHAR_ARG2 (&qz_job, 1),
432  F77_CONST_CHAR_ARG2 (&comp_q, 1),
433  F77_CONST_CHAR_ARG2 (&comp_z, 1),
434  nn, ilo, ihi,
437  F77_DBLE_CMPLX_ARG (xalpha.fortran_vec ()),
438  F77_DBLE_CMPLX_ARG (xbeta.fortran_vec ()),
439  F77_DBLE_CMPLX_ARG (CQ.fortran_vec ()), nn,
440  F77_DBLE_CMPLX_ARG (CZ.fortran_vec ()), nn,
441  F77_DBLE_CMPLX_ARG (cwork.fortran_vec ()), nn,
442  rwork.fortran_vec (), info
443  F77_CHAR_ARG_LEN (1)
444  F77_CHAR_ARG_LEN (1)
445  F77_CHAR_ARG_LEN (1)));
446 
447  if (comp_q == 'V')
448  {
449  // Left eigenvector.
450  F77_XFCN (zggbak, ZGGBAK,
451  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
452  F77_CONST_CHAR_ARG2 ("L", 1),
453  nn, ilo, ihi, lscale.data (), rscale.data (),
454  nn, F77_DBLE_CMPLX_ARG (CQ.fortran_vec ()), nn, info
455  F77_CHAR_ARG_LEN (1)
456  F77_CHAR_ARG_LEN (1)));
457  }
458 
459  if (comp_z == 'V')
460  {
461  // Right eigenvector.
462  F77_XFCN (zggbak, ZGGBAK,
463  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
464  F77_CONST_CHAR_ARG2 ("R", 1),
465  nn, ilo, ihi, lscale.data (), rscale.data (),
466  nn, F77_DBLE_CMPLX_ARG (CZ.fortran_vec ()), nn, info
467  F77_CHAR_ARG_LEN (1)
468  F77_CHAR_ARG_LEN (1)));
469  }
470 
471  }
472  else
473  {
474 #if defined (DEBUG)
475  std::cout << "qz: performing qr decomposition of bb" << std::endl;
476 #endif
477 
478  // Compute the QR factorization of bb.
479  octave::math::qr<Matrix> bqr (bb);
480 
481 #if defined (DEBUG)
482  std::cout << "qz: qr (bb) done; now performing qz decomposition"
483  << std::endl;
484 #endif
485 
486  bb = bqr.R ();
487 
488 #if defined (DEBUG)
489  std::cout << "qz: extracted bb" << std::endl;
490 #endif
491 
492  aa = (bqr.Q ()).transpose () * aa;
493 
494 #if defined (DEBUG)
495  std::cout << "qz: updated aa " << std::endl;
496  std::cout << "bqr.Q () = " << std::endl << bqr.Q () << std::endl;
497 
498  if (comp_q == 'V')
499  std::cout << "QQ =" << QQ << std::endl;
500 #endif
501 
502  if (comp_q == 'V')
503  QQ = QQ * bqr.Q ();
504 
505 #if defined (DEBUG)
506  std::cout << "qz: precursors done..." << std::endl;
507 #endif
508 
509 #if defined (DEBUG)
510  std::cout << "qz: comp_q = " << comp_q << ", comp_z = " << comp_z
511  << std::endl;
512 #endif
513 
514  // Reduce to generalized Hessenberg form.
515  F77_XFCN (dgghrd, DGGHRD,
516  (F77_CONST_CHAR_ARG2 (&comp_q, 1),
517  F77_CONST_CHAR_ARG2 (&comp_z, 1),
518  nn, ilo, ihi, aa.fortran_vec (),
519  nn, bb.fortran_vec (), nn, QQ.fortran_vec (), nn,
520  ZZ.fortran_vec (), nn, info
521  F77_CHAR_ARG_LEN (1)
522  F77_CHAR_ARG_LEN (1)));
523 
524  // Check if just computing generalized eigenvalues,
525  // or if we're actually computing the decomposition.
526 
527  // Reduce to generalized Schur form.
528  F77_XFCN (dhgeqz, DHGEQZ,
529  (F77_CONST_CHAR_ARG2 (&qz_job, 1),
530  F77_CONST_CHAR_ARG2 (&comp_q, 1),
531  F77_CONST_CHAR_ARG2 (&comp_z, 1),
532  nn, ilo, ihi, aa.fortran_vec (), nn, bb.fortran_vec (),
533  nn, alphar.fortran_vec (), alphai.fortran_vec (),
534  betar.fortran_vec (), QQ.fortran_vec (), nn,
535  ZZ.fortran_vec (), nn, work.fortran_vec (), nn, info
536  F77_CHAR_ARG_LEN (1)
537  F77_CHAR_ARG_LEN (1)
538  F77_CHAR_ARG_LEN (1)));
539 
540  if (comp_q == 'V')
541  {
542  F77_XFCN (dggbak, DGGBAK,
543  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
544  F77_CONST_CHAR_ARG2 ("L", 1),
545  nn, ilo, ihi, lscale.data (), rscale.data (),
546  nn, QQ.fortran_vec (), nn, info
547  F77_CHAR_ARG_LEN (1)
548  F77_CHAR_ARG_LEN (1)));
549 
550 #if defined (DEBUG)
551  if (comp_q == 'V')
552  std::cout << "qz: balancing done; QQ=" << std::endl
553  << QQ << std::endl;
554 #endif
555  }
556 
557  // then right
558  if (comp_z == 'V')
559  {
560  F77_XFCN (dggbak, DGGBAK,
561  (F77_CONST_CHAR_ARG2 (&bal_job, 1),
562  F77_CONST_CHAR_ARG2 ("R", 1),
563  nn, ilo, ihi, lscale.data (), rscale.data (),
564  nn, ZZ.fortran_vec (), nn, info
565  F77_CHAR_ARG_LEN (1)
566  F77_CHAR_ARG_LEN (1)));
567 
568 #if defined (DEBUG)
569  if (comp_z == 'V')
570  std::cout << "qz: balancing done; ZZ=" << std::endl
571  << ZZ << std::endl;
572 #endif
573  }
574 
575  }
576 
577  // Order the QZ decomposition?
578  if (ord_job != 'N')
579  {
580  if (complex_case)
581  // Probably not needed, but better be safe.
582  error ("qz: cannot re-order complex QZ decomposition");
583 
584 #if defined (DEBUG_SORT)
585  std::cout << "qz: ordering eigenvalues: ord_job = "
586  << ord_job << std::endl;
587 #endif
588 
589  Array<F77_LOGICAL> select (dim_vector (nn, 1));
590 
591  for (int j = 0; j < nn; j++)
592  {
593  switch (ord_job)
594  {
595  case 'S':
596  select(j) = alphar(j)*alphar(j) + alphai(j)*alphai(j) < betar(j)*betar(j);
597  break;
598 
599  case 'B':
600  select(j) = alphar(j)*alphar(j) + alphai(j)*alphai(j) >= betar(j)*betar(j);
601  break;
602 
603  case '+':
604  select(j) = alphar(j) * betar(j) >= 0;
605  break;
606 
607  case '-':
608  select(j) = alphar(j) * betar(j) < 0;
609  break;
610 
611  default:
612  // Invalid order option
613  // (should never happen since options were checked at the top).
614  panic_impossible ();
615  break;
616  }
617  }
618 
619  F77_LOGICAL wantq = 0, wantz = (comp_z == 'V');
620  F77_INT ijob = 0, mm, lrwork3 = 4*nn+16, liwork = nn;
621  F77_DBLE pl, pr;
622  RowVector rwork3(lrwork3);
623  Array<F77_INT> iwork (dim_vector (liwork, 1));
624 
625  F77_XFCN (dtgsen, DTGSEN,
626  (ijob, wantq, wantz,
627  select.fortran_vec (), nn,
628  aa.fortran_vec (), nn,
629  bb.fortran_vec (), nn,
630  alphar.fortran_vec (),
631  alphai.fortran_vec (),
632  betar.fortran_vec (),
633  nullptr, nn,
634  ZZ.fortran_vec (), nn,
635  mm,
636  pl, pr,
637  nullptr,
638  rwork3.fortran_vec (), lrwork3,
639  iwork.fortran_vec (), liwork,
640  info));
641 
642 #if defined (DEBUG_SORT)
643  std::cout << "qz: back from dtgsen: aa=" << std::endl;
644  octave_print_internal (std::cout, aa, 0);
645  std::cout << std::endl << "bb=" << std::endl;
646  octave_print_internal (std::cout, bb, 0);
647 
648  if (comp_z == 'V')
649  {
650  std::cout << std::endl << "ZZ=" << std::endl;
651  octave_print_internal (std::cout, ZZ, 0);
652  }
653  std::cout << std::endl << "qz: info=" << info << std::endl;
654  std::cout << "alphar = " << std::endl;
655  octave_print_internal (std::cout, (Matrix) alphar, 0);
656  std::cout << std::endl << "alphai = " << std::endl;
657  octave_print_internal (std::cout, (Matrix) alphai, 0);
658  std::cout << std::endl << "beta = " << std::endl;
659  octave_print_internal (std::cout, (Matrix) betar, 0);
660  std::cout << std::endl;
661 #endif
662  }
663 
664  // Compute the generalized eigenvalues as well?
666 
667  if (nargout < 2 || nargout == 7 || (nargin == 3 && nargout == 4))
668  {
669  if (complex_case)
670  {
672 
673  for (F77_INT i = 0; i < nn; i++)
674  tmp(i) = xalpha(i) / xbeta(i);
675 
676  gev = tmp;
677  }
678  else
679  {
680 #if defined (DEBUG)
681  std::cout << "qz: computing generalized eigenvalues" << std::endl;
682 #endif
683 
684  // Return finite generalized eigenvalues.
686 
687  F77_INT cnt = 0;
688  for (F77_INT i = 0; i < nn; i++)
689  if (betar(i) != 0)
690  tmp(cnt++) = Complex (alphar(i), alphai(i)) / betar(i);
691 
692  tmp.resize (cnt); // Trim vector to number of return values
693 
694  gev = tmp;
695  }
696  }
697 
698  // Right, left eigenvector matrices.
699  if (nargout >= 5)
700  {
701  // Which side to compute?
702  char side = (nargout == 5 ? 'R' : 'B');
703  // Compute all of them and backtransform
704  char howmany = 'B';
705  // Dummy pointer; select is not used.
706  F77_INT *select = nullptr;
707 
708  if (complex_case)
709  {
710  CVL = CQ;
711  CVR = CZ;
712  ComplexRowVector cwork2 (2 * nn);
713  RowVector rwork2 (8 * nn);
714  F77_INT m;
715 
716  F77_XFCN (ztgevc, ZTGEVC,
717  (F77_CONST_CHAR_ARG2 (&side, 1),
718  F77_CONST_CHAR_ARG2 (&howmany, 1),
719  select, nn, F77_DBLE_CMPLX_ARG (caa.fortran_vec ()), nn,
721  nn, F77_DBLE_CMPLX_ARG (CVL.fortran_vec ()), nn,
722  F77_DBLE_CMPLX_ARG (CVR.fortran_vec ()), nn, nn,
723  m, F77_DBLE_CMPLX_ARG (cwork2.fortran_vec ()),
724  rwork2.fortran_vec (), info
725  F77_CHAR_ARG_LEN (1)
726  F77_CHAR_ARG_LEN (1)));
727  }
728  else
729  {
730 #if defined (DEBUG)
731  std::cout << "qz: computing generalized eigenvectors" << std::endl;
732 #endif
733 
734  VL = QQ;
735  VR = ZZ;
736  F77_INT m;
737 
738  F77_XFCN (dtgevc, DTGEVC,
739  (F77_CONST_CHAR_ARG2 (&side, 1),
740  F77_CONST_CHAR_ARG2 (&howmany, 1),
741  select, nn, aa.fortran_vec (), nn, bb.fortran_vec (),
742  nn, VL.fortran_vec (), nn, VR.fortran_vec (), nn, nn,
743  m, work.fortran_vec (), info
744  F77_CHAR_ARG_LEN (1)
745  F77_CHAR_ARG_LEN (1)));
746 
747  // Now construct the complex form of VV, WW.
748  F77_INT j = 0;
749 
750  while (j < nn)
751  {
752  octave_quit ();
753 
754  // See if real or complex eigenvalue.
755 
756  // Column increment; assume complex eigenvalue.
757  int cinc = 2;
758 
759  if (j == (nn-1))
760  // Single column.
761  cinc = 1;
762  else if (aa(j+1,j) == 0)
763  cinc = 1;
764 
765  // Now copy the eigenvector (s) to CVR, CVL.
766  if (cinc == 1)
767  {
768  for (F77_INT i = 0; i < nn; i++)
769  CVR(i,j) = VR(i,j);
770 
771  if (side == 'B')
772  for (F77_INT i = 0; i < nn; i++)
773  CVL(i,j) = VL(i,j);
774  }
775  else
776  {
777  // Double column; complex vector.
778 
779  for (F77_INT i = 0; i < nn; i++)
780  {
781  CVR(i,j) = Complex (VR(i,j), VR(i,j+1));
782  CVR(i,j+1) = Complex (VR(i,j), -VR(i,j+1));
783  }
784 
785  if (side == 'B')
786  for (F77_INT i = 0; i < nn; i++)
787  {
788  CVL(i,j) = Complex (VL(i,j), VL(i,j+1));
789  CVL(i,j+1) = Complex (VL(i,j), -VL(i,j+1));
790  }
791  }
792 
793  // Advance to next eigenvectors (if any).
794  j += cinc;
795  }
796  }
797  }
798 
800 
801  switch (nargout)
802  {
803  case 7:
804  retval(6) = gev;
805  OCTAVE_FALLTHROUGH;
806 
807  case 6:
808  // Return left eigenvectors.
809  retval(5) = CVL;
810  OCTAVE_FALLTHROUGH;
811 
812  case 5:
813  // Return right eigenvectors.
814  retval(4) = CVR;
815  OCTAVE_FALLTHROUGH;
816 
817  case 4:
818  if (nargin == 3)
819  {
820 #if defined (DEBUG)
821  std::cout << "qz: sort: retval(3) = gev = " << std::endl;
822  octave_print_internal (std::cout, ComplexMatrix (gev));
823  std::cout << std::endl;
824 #endif
825  retval(3) = gev;
826  }
827  else
828  {
829  if (complex_case)
830  retval(3) = CZ;
831  else
832  retval(3) = ZZ;
833  }
834  OCTAVE_FALLTHROUGH;
835 
836  case 3:
837  if (nargin == 3)
838  {
839  if (complex_case)
840  retval(2) = CZ;
841  else
842  retval(2) = ZZ;
843  }
844  else
845  {
846  if (complex_case)
847  retval(2) = CQ.hermitian ();
848  else
849  retval(2) = QQ.transpose ();
850  }
851  OCTAVE_FALLTHROUGH;
852 
853  case 2:
854  {
855  if (complex_case)
856  {
857 #if defined (DEBUG)
858  std::cout << "qz: retval(1) = cbb = " << std::endl;
859  octave_print_internal (std::cout, cbb, 0);
860  std::cout << std::endl << "qz: retval(0) = caa = " <<std::endl;
861  octave_print_internal (std::cout, caa, 0);
862  std::cout << std::endl;
863 #endif
864  retval(1) = cbb;
865  retval(0) = caa;
866  }
867  else
868  {
869 #if defined (DEBUG)
870  std::cout << "qz: retval(1) = bb = " << std::endl;
871  octave_print_internal (std::cout, bb, 0);
872  std::cout << std::endl << "qz: retval(0) = aa = " <<std::endl;
873  octave_print_internal (std::cout, aa, 0);
874  std::cout << std::endl;
875 #endif
876  retval(1) = bb;
877  retval(0) = aa;
878  }
879  }
880  break;
881 
882  case 1:
883  case 0:
884 #if defined (DEBUG)
885  std::cout << "qz: retval(0) = gev = " << gev << std::endl;
886 #endif
887  retval(0) = gev;
888  break;
889 
890  default:
891  error ("qz: too many return arguments");
892  break;
893  }
894 
895 #if defined (DEBUG)
896  std::cout << "qz: exiting (at long last)" << std::endl;
897 #endif
898 
899  return retval;
900 }
901 
902 /*
903 %!shared a, b, c
904 %! a = [1 2; 0 3];
905 %! b = [1 0; 0 0];
906 %! c = [0 1; 0 0];
907 %!assert (qz (a,b), 1)
908 %!assert (isempty (qz (a,c)))
909 
910 ## Example 7.7.3 in Golub & Van Loan
911 %!test
912 %! a = [ 10 1 2;
913 %! 1 2 -1;
914 %! 1 1 2];
915 %! b = reshape (1:9,3,3);
916 %! [aa, bb, q, z, v, w, lambda] = qz (a, b);
917 %! sz = length (lambda);
918 %! observed = (b * v * diag ([lambda;0])) (:, 1:sz);
919 %! assert ((a*v)(:, 1:sz), observed, norm (observed) * 1e-14);
920 %! observed = (diag ([lambda;0]) * w' * b) (1:sz, :);
921 %! assert ((w'*a)(1:sz, :) , observed, norm (observed) * 1e-13);
922 %! assert (q * a * z, aa, norm (aa) * 1e-14);
923 %! assert (q * b * z, bb, norm (bb) * 1e-14);
924 
925 %!test
926 %! A = [0, 0, -1, 0; 1, 0, 0, 0; -1, 0, -2, -1; 0, -1, 1, 0];
927 %! B = [0, 0, 0, 0; 0, 1, 0, 0; 0, 0, 1, 0; 0, 0, 0, 1];
928 %! [AA, BB, Q, Z1] = qz (A, B);
929 %! [AA, BB, Z2] = qz (A, B, "-");
930 %! assert (Z1, Z2);
931 
932 %!test
933 %! A = [ -1.03428 0.24929 0.43205 -0.12860;
934 %! 1.16228 0.27870 2.12954 0.69250;
935 %! -0.51524 -0.34939 -0.77820 2.13721;
936 %! -1.32941 2.11870 0.72005 1.00835 ];
937 %! B = [ 1.407302 -0.632956 -0.360628 0.068534;
938 %! 0.149898 0.298248 0.991777 0.023652;
939 %! 0.169281 -0.405205 -1.775834 1.511730;
940 %! 0.717770 1.291390 -1.766607 -0.531352 ];
941 %! [AA, BB, Z, lambda] = qz (A, B, "+");
942 %! assert (all (real (lambda(1:3)) >= 0))
943 %! assert (real (lambda(4) < 0))
944 %! [AA, BB, Z, lambda] = qz (A, B, "-");
945 %! assert (real (lambda(1) < 0))
946 %! assert (all (real (lambda(2:4)) >= 0))
947 %! [AA, BB, Z, lambda] = qz (A, B, "B");
948 %! assert (all (abs (lambda(1:3)) >= 1))
949 %! assert (abs (lambda(4) < 1))
950 %! [AA, BB, Z, lambda] = qz (A, B, "S");
951 %! assert (abs (lambda(1) < 1))
952 %! assert (all (abs (lambda(2:4)) >= 1))
953 */
void warn_empty_arg(const char *name)
Definition: errwarn.cc:332
octave_idx_type rows(void) const
Definition: Array.h:404
F77_RET_T const F77_INT const F77_INT const F77_INT F77_DBLE_CMPLX const F77_INT F77_DBLE_CMPLX const F77_INT F77_DBLE_CMPLX F77_DBLE_CMPLX F77_DBLE_CMPLX const F77_INT F77_DBLE_CMPLX * CZ
subroutine xdlamch(cmach, retval)
Definition: xdlamch.f:2
bool isempty(void) const
Definition: Array.h:565
OCTINTERP_API void print_usage(void)
Definition: defun.cc:54
#define F77_DBLE_CMPLX_ARG(x)
Definition: f77-fcn.h:315
const T * fortran_vec(void) const
Definition: Array.h:584
static void transpose(octave_idx_type N, const octave_idx_type *ridx, const octave_idx_type *cidx, octave_idx_type *ridx2, octave_idx_type *cidx2)
Definition: symrcm.cc:386
#define DEFUN(name, args_name, nargout_name, doc)
Macro to define a builtin function.
Definition: defun.h:53
void error(const char *fmt,...)
Definition: error.cc:578
static F77_INT nn
Definition: DASPK.cc:62
octave_idx_type b_nr
Definition: sylvester.cc:76
octave_idx_type columns(void) const
Definition: Array.h:413
void err_square_matrix_required(const char *fcn, const char *name)
Definition: errwarn.cc:118
#define F77_XFCN(f, F, args)
Definition: f77-fcn.h:41
octave_f77_int_type F77_LOGICAL
Definition: f77-fcn.h:307
void octave_print_internal(std::ostream &os, const float_display_format &fmt, bool d, bool pr_as_read_syntax)
Definition: pr-output.cc:1780
OCTAVE_EXPORT octave_value_list return the number of command line arguments passed to Octave If called with the optional argument the function xample nargout(@histc)
Definition: ov-usr-fcn.cc:997
F77_RET_T F77_INT const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE * VL
double tmp
Definition: data.cc:6252
void err_nonconformant(const char *op, octave_idx_type op1_len, octave_idx_type op2_len)
octave_value retval
Definition: data.cc:6246
#define panic_impossible()
Definition: error.h:40
T Q(void) const
Definition: qr.h:78
Definition: dMatrix.h:36
F77_RET_T F77_FUNC(xerbla, XERBLA)
Definition: xerbla.c:57
F77_RET_T F77_INT const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE * VR
N Dimensional Array with copy-on-write semantics.
Definition: Array.h:125
octave_f77_int_type F77_INT
Definition: f77-fcn.h:305
octave_idx_type b_nc
Definition: sylvester.cc:77
args.length() nargin
Definition: file-io.cc:589
for i
Definition: data.cc:5264
std::complex< double > Complex
Definition: oct-cmplx.h:31
T R(void) const
Definition: qr.h:80
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:87
F77_RET_T const F77_INT const F77_INT const F77_INT F77_DBLE_CMPLX const F77_INT F77_DBLE_CMPLX const F77_INT F77_DBLE_CMPLX F77_DBLE_CMPLX F77_DBLE_CMPLX * CQ
If this string is the system will ring the terminal sometimes it is useful to be able to print the original representation of the string
Definition: utils.cc:888
ComplexMatrix hermitian(void) const
Definition: CMatrix.h:166