GNU Octave 11.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
qz.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1998-2026 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26// Generalized eigenvalue balancing via LAPACK
27
28// Originally written by A. S. Hodel <scotte@eng.auburn.edu>, but is
29// substantially different with the change to use LAPACK.
30
31#undef OCTAVE_QZ_DEBUG
32
33#if defined (HAVE_CONFIG_H)
34# include "config.h"
35#endif
36
37#include <cctype>
38#include <cmath>
39
40#include "f77-fcn.h"
41#include "lapack-proto.h"
42#include "qr.h"
43#include "quit.h"
44
45#include "defun.h"
46#include "error.h"
47#include "errwarn.h"
48#include "ovl.h"
49
50#if defined (OCTAVE_QZ_DEBUG)
51# include "pager.h"
52# include "pr-output.h"
53#endif
54
56
57DEFUN (qz, args, nargout,
58 doc: /* -*- texinfo -*-
59@deftypefn {} {[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}] =} qz (@var{A}, @var{B})
60@deftypefnx {} {[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}] =} qz (@var{A}, @var{B}, @var{opt})
61Compute the QZ@tie{}decomposition of a generalized eigenvalue problem.
62
63The generalized eigenvalue problem is defined as
64@tex
65$$A x = \lambda B x$$
66@end tex
67@ifnottex
68
69@math{A x = @var{lambda} B x}
70
71@end ifnottex
72
73There are two calling forms of the function:
74
75@enumerate
76@item @code{[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}] = qz (@var{A}, @var{B})}
77
78Compute the complex QZ@tie{}decomposition, generalized eigenvectors, and
79generalized eigenvalues.
80@tex
81$$ AA = Q \cdot A \cdot Z, BB = Q \cdot B \cdot Z $$
82$$ A \cdot V \cdot {\rm diag}(BB) = B \cdot V \cdot {\rm diag}(AA) $$
83$$ {\rm diag}(BB) \cdot W^T \cdot A = {\rm diag}(AA) \cdot W^T \cdot B $$
84@end tex
85@ifnottex
86
87@example
88@group
89
90@var{AA} = @var{Q} * @var{A} * @var{Z}, @var{BB} = @var{Q} * @var{B} * @var{Z}
91@var{A} * @var{V} * diag (diag (@var{BB})) = @var{B} * @var{V} * diag (diag (@var{AA}))
92diag (diag (@var{BB})) * @var{W}' * @var{A} = diag (diag (@var{AA})) * @var{W}' * @var{B}
93
94@end group
95@end example
96
97@end ifnottex
98with @var{AA} and @var{BB} upper triangular, and @var{Q} and @var{Z}
99unitary. The matrices @var{V} and @var{W} respectively contain the right
100and left generalized eigenvectors.
101
102@item @code{[@var{AA}, @var{BB}, @var{Q}, @var{Z}, @var{V}, @var{W}] = qz (@var{A}, @var{B}, @var{opt})}
103
104The @var{opt} argument must be equal to either @qcode{"real"} or
105@qcode{"complex"}. If it is equal to @qcode{"complex"}, then this
106calling form is equivalent to the first one with only two input
107arguments.
108
109If @var{opt} is equal to @qcode{"real"}, then the real QZ decomposition
110is computed. In particular, @var{AA} is only guaranteed to be
111quasi-upper triangular with 1-by-1 and 2-by-2 blocks on the diagonal,
112and @var{Q} and @var{Z} are orthogonal. The identities mentioned above
113for right and left generalized eigenvectors are only verified if
114@var{AA} is upper triangular (i.e., when all the generalized eigenvalues
115are real, in which case the real and complex QZ coincide).
116
117@end enumerate
118
119Note: @code{qz} performs permutation balancing, but not scaling
120(@pxref{XREFbalance,,@code{balance}}), which may lead to less accurate
121results than @code{eig}. The order of output arguments was selected for
122compatibility with @sc{matlab}.
123
124For eigenvalues, use the @code{ordeig} function to compute them based on
125the resulting @var{AA} and @var{BB} matrices.
126
127@seealso{eig, gsvd, balance, chol, hess, lu, qr, qzhess, schur, ordeig}
128@end deftypefn */)
129{
130 int nargin = args.length ();
131
132#if defined (OCTAVE_QZ_DEBUG)
133 octave_stdout << "qz: nargin = " << nargin
134 << ", nargout = " << nargout << std::endl;
135#endif
136
137 if (nargin < 2 || nargin > 3 || nargout > 6)
138 print_usage ();
139
140 bool complex_case = true;
141
142 if (nargin == 3)
143 {
144 std::string opt = args(2).xstring_value ("qz: OPT must be a string");
145
146 if (opt == "real")
147 {
148 if (args(0).iscomplex () || args(1).iscomplex ())
149 {
150 warning ("qz: ignoring 'real' option with complex matrices");
151 complex_case = true;
152 }
153 else
154 complex_case = false;
155 }
156 else if (opt == "complex")
157 complex_case = true;
158 else
159 error ("qz: OPT must be 'real' or 'complex'");
160 }
161
162#if defined (OCTAVE_QZ_DEBUG)
163 octave_stdout << "qz: check matrix A" << std::endl;
164#endif
165
166 // Matrix A: check dimensions.
167 F77_INT nn = to_f77_int (args(0).rows ());
168 F77_INT nc = to_f77_int (args(0).columns ());
169
170#if defined (OCTAVE_QZ_DEBUG)
171 octave_stdout << "Matrix A dimensions: (" << nn << ',' << nc << ')' << std::endl;
172#endif
173
174 if (nc != nn)
175 err_square_matrix_required ("qz", "A");
176
177 // Matrix A: get value.
178 Matrix aa;
179 ComplexMatrix caa;
180
181 if (args(0).iscomplex ())
182 caa = args(0).complex_matrix_value ();
183 else
184 aa = args(0).matrix_value ();
185
186#if defined (OCTAVE_QZ_DEBUG)
187 octave_stdout << "qz: check matrix B" << std::endl;
188#endif
189
190 // Extract argument 2 (bb, or cbb if complex).
191 F77_INT b_nr = to_f77_int (args(1).rows ());
192 F77_INT b_nc = to_f77_int (args(1).columns ());
193
194 if (nn != b_nc || nn != b_nr)
196
197 Matrix bb;
198 ComplexMatrix cbb;
199
200 if (args(1).iscomplex ())
201 cbb = args(1).complex_matrix_value ();
202 else
203 bb = args(1).matrix_value ();
204
205 // First, declare variables used in both the real and complex cases.
206 // FIXME: There are a lot of excess variables declared.
207 // Probably a better way to handle this.
208 Matrix QQ (nn, nn), ZZ (nn, nn), VR (nn, nn), VL (nn, nn);
209 RowVector alphar (nn), alphai (nn), betar (nn);
210 ComplexRowVector xalpha (nn), xbeta (nn);
211 ComplexMatrix CQ (nn, nn), CZ (nn, nn), CVR (nn, nn), CVL (nn, nn);
212 F77_INT ilo, ihi, info;
213 char comp_q = (nargout >= 3 ? 'V' : 'N');
214 char comp_z = (nargout >= 4 ? 'V' : 'N');
215
216 // Initialize Q, Z to identity matrix if either is needed
217 if (comp_q == 'V' || comp_z == 'V')
218 {
219 double *QQptr = QQ.rwdata ();
220 double *ZZptr = ZZ.rwdata ();
221 std::fill_n (QQptr, QQ.numel (), 0.0);
222 std::fill_n (ZZptr, ZZ.numel (), 0.0);
223 for (F77_INT i = 0; i < nn; i++)
224 {
225 QQ(i, i) = 1.0;
226 ZZ(i, i) = 1.0;
227 }
228 }
229
230 // Always perform permutation balancing.
231 const char bal_job = 'P';
232 RowVector lscale (nn), rscale (nn), work (6 * nn), rwork (nn);
233
234 if (complex_case)
235 {
236 if (args(0).isreal ())
237 caa = ComplexMatrix (aa);
238
239 if (args(1).isreal ())
240 cbb = ComplexMatrix (bb);
241
242 if (comp_q == 'V')
243 CQ = ComplexMatrix (QQ);
244
245#if defined (OCTAVE_QZ_DEBUG)
246 if (comp_q == 'V')
247 octave_stdout << "qz: performing balancing; CQ =\n" << CQ << std::endl;
248#endif
249
250 if (comp_z == 'V')
251 CZ = ComplexMatrix (ZZ);
252
253 F77_XFCN (zggbal, ZGGBAL,
254 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
255 nn, F77_DBLE_CMPLX_ARG (caa.rwdata ()), nn,
256 F77_DBLE_CMPLX_ARG (cbb.rwdata ()),
257 nn, ilo, ihi, lscale.rwdata (),
258 rscale.rwdata (), work.rwdata (), info
259 F77_CHAR_ARG_LEN (1)));
260 }
261 else
262 {
263#if defined (OCTAVE_QZ_DEBUG)
264 if (comp_q == 'V')
265 octave_stdout << "qz: performing balancing; QQ =\n" << QQ << std::endl;
266#endif
267
268 F77_XFCN (dggbal, DGGBAL,
269 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
270 nn, aa.rwdata (), nn, bb.rwdata (),
271 nn, ilo, ihi, lscale.rwdata (),
272 rscale.rwdata (), work.rwdata (), info
273 F77_CHAR_ARG_LEN (1)));
274 }
275
276 // Only permutation balance above is done. Skip scaling balance.
277
278 char qz_job = (nargout < 2 ? 'E' : 'S');
279
280 if (complex_case)
281 {
282 // Complex case.
283
284 // The QR decomposition of cbb.
285 math::qr<ComplexMatrix> cbqr (cbb);
286 // The R matrix of QR decomposition for cbb.
287 cbb = cbqr.R ();
288 // (Q*)caa for following work.
289 caa = (cbqr.Q ().hermitian ()) * caa;
290 CQ = CQ * cbqr.Q ();
291
292 F77_XFCN (zgghrd, ZGGHRD,
293 (F77_CONST_CHAR_ARG2 (&comp_q, 1),
294 F77_CONST_CHAR_ARG2 (&comp_z, 1),
295 nn, ilo, ihi, F77_DBLE_CMPLX_ARG (caa.rwdata ()),
296 nn, F77_DBLE_CMPLX_ARG (cbb.rwdata ()), nn,
297 F77_DBLE_CMPLX_ARG (CQ.rwdata ()), nn,
298 F77_DBLE_CMPLX_ARG (CZ.rwdata ()), nn, info
299 F77_CHAR_ARG_LEN (1)
300 F77_CHAR_ARG_LEN (1)));
301
302 ComplexRowVector cwork (nn);
303
304 F77_XFCN (zhgeqz, ZHGEQZ,
305 (F77_CONST_CHAR_ARG2 (&qz_job, 1),
306 F77_CONST_CHAR_ARG2 (&comp_q, 1),
307 F77_CONST_CHAR_ARG2 (&comp_z, 1),
308 nn, ilo, ihi,
309 F77_DBLE_CMPLX_ARG (caa.rwdata ()), nn,
310 F77_DBLE_CMPLX_ARG (cbb.rwdata ()), nn,
311 F77_DBLE_CMPLX_ARG (xalpha.rwdata ()),
312 F77_DBLE_CMPLX_ARG (xbeta.rwdata ()),
313 F77_DBLE_CMPLX_ARG (CQ.rwdata ()), nn,
314 F77_DBLE_CMPLX_ARG (CZ.rwdata ()), nn,
315 F77_DBLE_CMPLX_ARG (cwork.rwdata ()), nn,
316 rwork.rwdata (), info
317 F77_CHAR_ARG_LEN (1)
318 F77_CHAR_ARG_LEN (1)
319 F77_CHAR_ARG_LEN (1)));
320
321 if (comp_q == 'V')
322 {
323 // Left eigenvector.
324 F77_XFCN (zggbak, ZGGBAK,
325 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
326 F77_CONST_CHAR_ARG2 ("L", 1),
327 nn, ilo, ihi, lscale.data (), rscale.data (),
328 nn, F77_DBLE_CMPLX_ARG (CQ.rwdata ()), nn, info
329 F77_CHAR_ARG_LEN (1)
330 F77_CHAR_ARG_LEN (1)));
331 }
332
333 if (comp_z == 'V')
334 {
335 // Right eigenvector.
336 F77_XFCN (zggbak, ZGGBAK,
337 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
338 F77_CONST_CHAR_ARG2 ("R", 1),
339 nn, ilo, ihi, lscale.data (), rscale.data (),
340 nn, F77_DBLE_CMPLX_ARG (CZ.rwdata ()), nn, info
341 F77_CHAR_ARG_LEN (1)
342 F77_CHAR_ARG_LEN (1)));
343 }
344
345 }
346 else
347 {
348#if defined (OCTAVE_QZ_DEBUG)
349 octave_stdout << "qz: performing qr decomposition of bb" << std::endl;
350#endif
351
352 // Compute the QR factorization of bb.
353 math::qr<Matrix> bqr (bb);
354
355#if defined (OCTAVE_QZ_DEBUG)
356 octave_stdout << "qz: qr (bb) done; now performing qz decomposition" << std::endl;
357#endif
358
359 bb = bqr.R ();
360
361#if defined (OCTAVE_QZ_DEBUG)
362 octave_stdout << "qz: extracted bb" << std::endl;
363#endif
364
365 aa = (bqr.Q ()).transpose () * aa;
366
367#if defined (OCTAVE_QZ_DEBUG)
368 octave_stdout << "qz: updated aa " << std::endl;
369 octave_stdout << "bqr.Q () =\n" << bqr.Q () << std::endl;
370
371 if (comp_q == 'V')
372 octave_stdout << "QQ =" << QQ << std::endl;
373#endif
374
375 if (comp_q == 'V')
376 QQ = QQ * bqr.Q ();
377
378#if defined (OCTAVE_QZ_DEBUG)
379 octave_stdout << "qz: precursors done..." << std::endl;
380 octave_stdout << "qz: comp_q = " << comp_q << ", comp_z = " << comp_z << std::endl;
381#endif
382
383 // Reduce to generalized Hessenberg form.
384 F77_XFCN (dgghrd, DGGHRD,
385 (F77_CONST_CHAR_ARG2 (&comp_q, 1),
386 F77_CONST_CHAR_ARG2 (&comp_z, 1),
387 nn, ilo, ihi, aa.rwdata (),
388 nn, bb.rwdata (), nn, QQ.rwdata (), nn,
389 ZZ.rwdata (), nn, info
390 F77_CHAR_ARG_LEN (1)
391 F77_CHAR_ARG_LEN (1)));
392
393 // Check if just computing generalized eigenvalues,
394 // or if we're actually computing the decomposition.
395
396 // Reduce to generalized Schur form.
397 F77_XFCN (dhgeqz, DHGEQZ,
398 (F77_CONST_CHAR_ARG2 (&qz_job, 1),
399 F77_CONST_CHAR_ARG2 (&comp_q, 1),
400 F77_CONST_CHAR_ARG2 (&comp_z, 1),
401 nn, ilo, ihi, aa.rwdata (), nn, bb.rwdata (),
402 nn, alphar.rwdata (), alphai.rwdata (),
403 betar.rwdata (), QQ.rwdata (), nn,
404 ZZ.rwdata (), nn, work.rwdata (), nn, info
405 F77_CHAR_ARG_LEN (1)
406 F77_CHAR_ARG_LEN (1)
407 F77_CHAR_ARG_LEN (1)));
408
409 if (comp_q == 'V')
410 {
411 F77_XFCN (dggbak, DGGBAK,
412 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
413 F77_CONST_CHAR_ARG2 ("L", 1),
414 nn, ilo, ihi, lscale.data (), rscale.data (),
415 nn, QQ.rwdata (), nn, info
416 F77_CHAR_ARG_LEN (1)
417 F77_CHAR_ARG_LEN (1)));
418
419#if defined (OCTAVE_QZ_DEBUG)
420 if (comp_q == 'V')
421 octave_stdout << "qz: balancing done; QQ=\n" << QQ << std::endl;
422#endif
423 }
424
425 // then right
426 if (comp_z == 'V')
427 {
428 F77_XFCN (dggbak, DGGBAK,
429 (F77_CONST_CHAR_ARG2 (&bal_job, 1),
430 F77_CONST_CHAR_ARG2 ("R", 1),
431 nn, ilo, ihi, lscale.data (), rscale.data (),
432 nn, ZZ.rwdata (), nn, info
433 F77_CHAR_ARG_LEN (1)
434 F77_CHAR_ARG_LEN (1)));
435
436#if defined (OCTAVE_QZ_DEBUG)
437 if (comp_z == 'V')
438 octave_stdout << "qz: balancing done; ZZ=\n" << ZZ << std::endl;
439#endif
440 }
441
442 }
443
444 // Right, left eigenvector matrices.
445 if (nargout >= 5)
446 {
447 // Which side to compute?
448 char side = (nargout == 5 ? 'R' : 'B');
449 // Compute all of them and backtransform
450 char howmany = 'B';
451 // Dummy pointer; select is not used.
452 F77_INT *select = nullptr;
453
454 if (complex_case)
455 {
456 CVL = CQ;
457 CVR = CZ;
458 ComplexRowVector cwork2 (2 * nn);
459 RowVector rwork2 (8 * nn);
460 F77_INT m;
461
462 F77_XFCN (ztgevc, ZTGEVC,
463 (F77_CONST_CHAR_ARG2 (&side, 1),
464 F77_CONST_CHAR_ARG2 (&howmany, 1),
465 select, nn, F77_DBLE_CMPLX_ARG (caa.rwdata ()), nn,
466 F77_DBLE_CMPLX_ARG (cbb.rwdata ()),
467 nn, F77_DBLE_CMPLX_ARG (CVL.rwdata ()), nn,
468 F77_DBLE_CMPLX_ARG (CVR.rwdata ()), nn, nn,
469 m, F77_DBLE_CMPLX_ARG (cwork2.rwdata ()),
470 rwork2.rwdata (), info
471 F77_CHAR_ARG_LEN (1)
472 F77_CHAR_ARG_LEN (1)));
473 }
474 else
475 {
476#if defined (OCTAVE_QZ_DEBUG)
477 octave_stdout << "qz: computing generalized eigenvectors" << std::endl;
478#endif
479
480 VL = QQ;
481 VR = ZZ;
482 F77_INT m;
483
484 F77_XFCN (dtgevc, DTGEVC,
485 (F77_CONST_CHAR_ARG2 (&side, 1),
486 F77_CONST_CHAR_ARG2 (&howmany, 1),
487 select, nn, aa.rwdata (), nn, bb.rwdata (),
488 nn, VL.rwdata (), nn, VR.rwdata (), nn, nn,
489 m, work.rwdata (), info
490 F77_CHAR_ARG_LEN (1)
491 F77_CHAR_ARG_LEN (1)));
492 }
493 }
494
495 octave_value_list retval (nargout);
496
497 switch (nargout)
498 {
499 case 6:
500 // Return left eigenvectors.
501 if (complex_case)
502 retval(5) = CVL;
503 else
504 retval(5) = VL;
505 OCTAVE_FALLTHROUGH;
506
507 case 5:
508 // Return right eigenvectors.
509 if (complex_case)
510 retval(4) = CVR;
511 else
512 retval(4) = VR;
513 OCTAVE_FALLTHROUGH;
514
515 case 4:
516 if (complex_case)
517 retval(3) = CZ;
518 else
519 retval(3) = ZZ;
520 OCTAVE_FALLTHROUGH;
521
522 case 3:
523 if (complex_case)
524 retval(2) = CQ.hermitian ();
525 else
526 retval(2) = QQ.transpose ();
527 OCTAVE_FALLTHROUGH;
528
529 case 2:
530 case 1:
531 case 0:
532 {
533 if (complex_case)
534 {
535#if defined (OCTAVE_QZ_DEBUG)
536 octave_stdout << "qz: retval(1) = cbb =\n";
538 octave_stdout << "\nqz: retval(0) = caa =\n";
540 octave_stdout << std::endl;
541#endif
542 retval(1) = cbb;
543 retval(0) = caa;
544 }
545 else
546 {
547#if defined (OCTAVE_QZ_DEBUG)
548 octave_stdout << "qz: retval(1) = bb =\n";
550 octave_stdout << "\nqz: retval(0) = aa =\n";
552 octave_stdout << std::endl;
553#endif
554 retval(1) = bb;
555 retval(0) = aa;
556 }
557 }
558 break;
559
560 }
561
562#if defined (OCTAVE_QZ_DEBUG)
563 octave_stdout << "qz: exiting (at long last)" << std::endl;
564#endif
565
566 return retval;
567}
568
569/*
570## Example 7.7.3 in Golub & Van Loan (3rd ed.; not present in 4th ed.)
571## The generalized eigenvalues are real, so the real and complex QZ coincide.
572%!test
573%! a = [ 10 1 2;
574%! 1 2 -1;
575%! 1 1 2 ];
576%! b = reshape (1:9, 3,3);
577%! [aa, bb, q, z, v, w] = qz (a, b);
578%! assert (isreal (aa) && isreal (bb) && isreal (q) && isreal (z) ...
579%! && isreal (v) && isreal (w));
580%! assert (istriu (aa) && istriu (bb));
581%! assert (q * q', eye(3), 1e-15);
582%! assert (z * z', eye(3), 1e-15);
583%! assert (q * a * z, aa, norm (aa) * 1e-14);
584%! assert (q * b * z, bb, norm (bb) * 1e-14);
585%! assert (a * v * diag (diag (bb)), b * v * diag (diag (aa)), 1e-13);
586%! assert (diag (diag (bb)) * w' * a, diag (diag (aa)) * w' * b, 1e-13);
587
588## An example with real matrices where some generalized eigenvalues are
589## complex, so the real and complex QZ differ.
590%!test
591%! A = [ -1.03428 0.24929 0.43205 -0.12860;
592%! 1.16228 0.27870 2.12954 0.69250;
593%! -0.51524 -0.34939 -0.77820 2.13721;
594%! -1.32941 2.11870 0.72005 1.00835 ];
595%! B = [ 1.407302 -0.632956 -0.360628 0.068534;
596%! 0.149898 0.298248 0.991777 0.023652;
597%! 0.169281 -0.405205 -1.775834 1.511730;
598%! 0.717770 1.291390 -1.766607 -0.531352 ];
599%! [CAA, CBB, CQ, CZ, CV, CW] = qz (A, B, 'complex');
600%! assert (iscomplex (CAA) && iscomplex (CBB) && iscomplex (CQ) ...
601%! && iscomplex (CZ) && iscomplex (CV) && iscomplex (CW));
602%! assert (istriu (CAA) && istriu (CBB));
603%! assert (CQ * CQ', eye(4), 1e-14);
604%! assert (CZ * CZ', eye(4), 1e-14);
605%! assert (CQ * A * CZ, CAA, norm (CAA) * 1e-14);
606%! assert (CQ * B * CZ, CBB, norm (CBB) * 1e-14);
607%! assert (A * CV * diag (diag (CBB)), B * CV * diag (diag (CAA)), 1e-13);
608%! assert (diag (diag (CBB)) * CW' * A, diag (diag (CAA)) * CW' * B, 1e-13);
609%! [AA, BB, Q, Z, V, W] = qz (A, B, 'real');
610%! assert (isreal (AA) && isreal (BB) && isreal (Q) && isreal (Z) ...
611%! && isreal (V) && isreal (W));
612%! assert (istriu (BB));
613%! assert (Q * Q', eye(4), 1e-14);
614%! assert (Z * Z', eye(4), 1e-14);
615%! assert (Q * A * Z, AA, norm (AA) * 1e-14);
616%! assert (Q * B * Z, BB, norm (BB) * 1e-14);
617
618## Test input validation
619%!error <Invalid call> qz ()
620%!error <Invalid call> qz (1)
621%!error <Invalid call> qz (1,2,3,4)
622%!error <Invalid call> [r1,r2,r3,r4,r5,r6,r7] = qz (1,2)
623%!error <OPT must be a string> qz (1,2, 3)
624%!error <OPT must be 'real' or 'complex'> qz (1,2, 'foobar')
625%!warning <ignoring 'real' option with complex matrices> qz (2i, 3, 'real');
626%!warning <ignoring 'real' option with complex matrices> qz (2, 3i, 'real');
627
628*/
629
630OCTAVE_END_NAMESPACE(octave)
T * rwdata()
Size of the specified dimension.
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
void print_usage()
Definition defun-int.h:72
#define DEFUN(name, args_name, nargout_name, doc)
Macro to define a builtin function.
Definition defun.h:56
void warning(const char *fmt,...)
Definition error.cc:1083
void error(const char *fmt,...)
Definition error.cc:1008
void err_square_matrix_required(const char *fcn, const char *name)
Definition errwarn.cc:122
void err_nonconformant()
Definition errwarn.cc:95
#define F77_DBLE_CMPLX_ARG(x)
Definition f77-fcn.h:316
#define F77_XFCN(f, F, args)
Definition f77-fcn.h:45
octave_f77_int_type F77_INT
Definition f77-fcn.h:306
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
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
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
F77_RET_T F77_INT const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE * VL
#define octave_stdout
Definition pager.h:301
void octave_print_internal(std::ostream &os, const float_display_format &fmt, bool d, bool pr_as_read_syntax)