Mercurial > octave
annotate scripts/linear-algebra/qzhess.m @ 30564:796f54d4ddbf stable
update Octave Project Developers copyright for the new year
In files that have the "Octave Project Developers" copyright notice,
update for 2021.
In all .txi and .texi files except gpl.txi and gpl.texi in the
doc/liboctave and doc/interpreter directories, change the copyright
to "Octave Project Developers", the same as used for other source
files. Update copyright notices for 2022 (not done since 2019). For
gpl.txi and gpl.texi, change the copyright notice to be "Free Software
Foundation, Inc." and leave the date at 2007 only because this file
only contains the text of the GPL, not anything created by the Octave
Project Developers.
Add Paul Thomas to contributors.in.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Tue, 28 Dec 2021 18:22:40 -0500 |
parents | 7854d5752dd2 |
children | 597f3ee61a48 |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 ######################################################################## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 ## |
30564
796f54d4ddbf
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
29359
diff
changeset
|
3 ## Copyright (C) 1993-2022 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26600
diff
changeset
|
4 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 ## See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 ## distribution or <https://octave.org/copyright/>. |
2313 | 7 ## |
8 ## This file is part of Octave. | |
9 ## | |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
10 ## Octave is free software: you can redistribute it and/or modify it |
2313 | 11 ## under the terms of the GNU General Public License as published by |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
12 ## the Free Software Foundation, either version 3 of the License, or |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
13 ## (at your option) any later version. |
2313 | 14 ## |
15 ## Octave is distributed in the hope that it will be useful, but | |
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
18 ## GNU General Public License for more details. |
2313 | 19 ## |
20 ## You should have received a copy of the GNU General Public License | |
7016 | 21 ## along with Octave; see the file COPYING. If not, see |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
22 ## <https://www.gnu.org/licenses/>. |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
23 ## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 ######################################################################## |
245 | 25 |
3372 | 26 ## -*- texinfo -*- |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20160
diff
changeset
|
27 ## @deftypefn {} {[@var{aa}, @var{bb}, @var{q}, @var{z}] =} qzhess (@var{A}, @var{B}) |
3372 | 28 ## Compute the Hessenberg-triangular decomposition of the matrix pencil |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
29 ## @code{(@var{A}, @var{B})}, returning |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
30 ## @code{@var{aa} = @var{q} * @var{A} * @var{z}}, |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
31 ## @code{@var{bb} = @var{q} * @var{B} * @var{z}}, with @var{q} and @var{z} |
20160
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
32 ## orthogonal. |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
33 ## |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
34 ## For example: |
3426 | 35 ## |
3372 | 36 ## @example |
37 ## @group | |
38 ## [aa, bb, q, z] = qzhess ([1, 2; 3, 4], [5, 6; 7, 8]) | |
26600
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
39 ## @result{} aa = |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
40 ## -3.02244 -4.41741 |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
41 ## 0.92998 0.69749 |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
42 ## @result{} bb = |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
43 ## -8.60233 -9.99730 |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
44 ## 0.00000 -0.23250 |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
45 ## @result{} q = |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
46 ## -0.58124 -0.81373 |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
47 ## -0.81373 0.58124 |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
48 ## @result{} z = |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
49 ## Diagonal Matrix |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
50 ## 1 0 |
f6730533820e
doc: clean up doc example blocks for accuracy and consistent formatting
Mike Miller <mtmiller@octave.org>
parents:
26376
diff
changeset
|
51 ## 0 1 |
3372 | 52 ## @end group |
53 ## @end example | |
3426 | 54 ## |
3372 | 55 ## The Hessenberg-triangular decomposition is the first step in |
19040
0850b5212619
doc: Add @nospell macro around proper names in documentation.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
56 ## @nospell{Moler and Stewart's} QZ@tie{}decomposition algorithm. |
3426 | 57 ## |
19040
0850b5212619
doc: Add @nospell macro around proper names in documentation.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
58 ## Algorithm taken from @nospell{Golub and Van Loan}, |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
59 ## @cite{Matrix Computations, 2nd edition}. |
19593
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
60 ## |
16920
53eaa83e4181
doc: Add seealso links between various factorization forms.
Rik <rik@octave.org>
parents:
14363
diff
changeset
|
61 ## @seealso{lu, chol, hess, qr, qz, schur, svd} |
3372 | 62 ## @end deftypefn |
29 | 63 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
64 function [aa, bb, q, z] = qzhess (A, B) |
72 | 65 |
54 | 66 if (nargin != 2) |
6046 | 67 print_usage (); |
54 | 68 endif |
69 | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
70 [na, ma] = size (A); |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
71 [nb, mb] = size (B); |
54 | 72 if (na != ma || na != nb || nb != mb) |
73 error ("qzhess: incompatible dimensions"); | |
74 endif | |
75 | |
2303 | 76 ## Reduce to hessenberg-triangular form. |
29 | 77 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
78 [q, bb] = qr (B); |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
79 aa = q' * A; |
54 | 80 q = q'; |
81 z = eye (na); | |
82 for j = 1:(na-2) | |
83 for i = na:-1:(j+2) | |
84 | |
2303 | 85 ## disp (["zero out aa(", num2str(i), ",", num2str(j), ")"]) |
54 | 86 |
87 rot = givens (aa (i-1, j), aa (i, j)); | |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
88 aa((i-1):i, :) = rot *aa((i-1):i, :); |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
89 bb((i-1):i, :) = rot *bb((i-1):i, :); |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
90 q((i-1):i, :) = rot * q((i-1):i, :); |
54 | 91 |
2303 | 92 ## disp (["now zero out bb(", num2str(i), ",", num2str(i-1), ")"]) |
54 | 93 |
94 rot = givens (bb (i, i), bb (i, i-1))'; | |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
95 bb(:, (i-1):i) = bb(:, (i-1):i) * rot'; |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
96 aa(:, (i-1):i) = aa(:, (i-1):i) * rot'; |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
97 z(:, (i-1):i) = z(:, (i-1):i) * rot'; |
54 | 98 |
99 endfor | |
29 | 100 endfor |
101 | |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
102 bb(2, 1) = 0.0; |
54 | 103 for i = 3:na |
104 bb (i, 1:(i-1)) = zeros (1, i-1); | |
105 aa (i, 1:(i-2)) = zeros (1, i-2); | |
106 endfor | |
107 | |
29 | 108 endfunction |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
109 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
110 |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
111 %!test |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
112 %! a = [1 2 1 3; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
113 %! 2 5 3 2; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
114 %! 5 5 1 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
115 %! 4 0 3 2]; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
116 %! b = [0 4 2 1; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
117 %! 2 3 1 1; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
118 %! 1 0 2 1; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
119 %! 2 5 3 2]; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
120 %! mask = [0 0 0 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
121 %! 0 0 0 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
122 %! 1 0 0 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
123 %! 1 1 0 0]; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
124 %! [aa, bb, q, z] = qzhess (a, b); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
125 %! assert (inv (q) - q', zeros (4), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
126 %! assert (inv (z) - z', zeros (4), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
127 %! assert (q * a * z, aa, 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
128 %! assert (aa .* mask, zeros (4), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
129 %! assert (q * b * z, bb, 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
130 %! assert (bb .* mask, zeros (4), 2e-8); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
131 |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
132 %!test |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
133 %! a = [1 2 3 4 5; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
134 %! 3 2 3 1 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
135 %! 4 3 2 1 1; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
136 %! 0 1 0 1 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
137 %! 3 2 1 0 5]; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
138 %! b = [5 0 4 0 1; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
139 %! 1 1 1 2 5; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
140 %! 0 3 2 1 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
141 %! 4 3 0 3 5; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
142 %! 2 1 2 1 3]; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
143 %! mask = [0 0 0 0 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
144 %! 0 0 0 0 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
145 %! 1 0 0 0 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
146 %! 1 1 0 0 0; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
147 %! 1 1 1 0 0]; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
148 %! [aa, bb, q, z] = qzhess (a, b); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
149 %! assert (inv (q) - q', zeros (5), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
150 %! assert (inv (z) - z', zeros (5), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
151 %! assert (q * a * z, aa, 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
152 %! assert (aa .* mask, zeros (5), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
153 %! assert (q * b * z, bb, 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
154 %! assert (bb .* mask, zeros (5), 2e-8); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
11593
diff
changeset
|
155 |
28896
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
27985
diff
changeset
|
156 ## Test input validation |
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
27985
diff
changeset
|
157 %!error <Invalid call> qzhess () |
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
27985
diff
changeset
|
158 %!error <Invalid call> qzhess (1) |