Mercurial > forge
comparison main/linear-algebra/inst/condeig.m @ 2719:4e0bad669780 octave-forge
Initial commit into CVS.
author | whyly |
---|---|
date | Tue, 17 Oct 2006 22:18:07 +0000 |
parents | |
children | a33bfad6a8d9 |
comparison
equal
deleted
inserted
replaced
2718:78f956ba9369 | 2719:4e0bad669780 |
---|---|
1 ## Copyright (C) 2006 Arno Onken | |
2 ## | |
3 ## This program is free software; you can redistribute it and/or modify | |
4 ## it under the terms of the GNU General Public License as published by | |
5 ## the Free Software Foundation; either version 2 of the License, or | |
6 ## (at your option) any later version. | |
7 ## | |
8 ## This program is distributed in the hope that it will be useful, | |
9 ## but WITHOUT ANY WARRANTY; without even the implied warranty of | |
10 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the | |
11 ## GNU General Public License for more details. | |
12 ## | |
13 ## You should have received a copy of the GNU General Public License | |
14 ## along with this program; if not, write to the Free Software | |
15 ## Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA | |
16 | |
17 ## -*- texinfo -*- | |
18 ## @deftypefn {Function File} {@var{c} =} condeig (@var{a}) | |
19 ## @deftypefnx {Function File} {[@var{v}, @var{lambda}, @var{c}] =} condeig (@var{a}) | |
20 ## Computes condition numbers for the eigenvalues of a matrix. The | |
21 ## condition numbers are the reciprocals of the cosines of the angles | |
22 ## between the left and right eigenvectors. | |
23 ## | |
24 ## Arguments are | |
25 ## | |
26 ## @itemize @bullet | |
27 ## @item | |
28 ## @var{a} must be a square numeric matrix. | |
29 ## @end itemize | |
30 ## | |
31 ## Return values are | |
32 ## | |
33 ## @itemize @bullet | |
34 ## @item | |
35 ## @var{c} is a vector of condition numbers for the eigenvalue of | |
36 ## @var{a}. | |
37 ## | |
38 ## @item | |
39 ## @var{v} is the matrix of right eigenvectors of @var{a}. The result is | |
40 ## the same as for @code{[v, lambda] = eig (a)}. | |
41 ## | |
42 ## @item | |
43 ## @var{lambda} is the diagonal matrix of eigenvalues of @var{a}. The | |
44 ## result is the same as for @code{[v, lambda] = eig (a)}. | |
45 ## @end itemize | |
46 ## | |
47 ## Example: | |
48 ## | |
49 ## @example | |
50 ## @group | |
51 ## a = [1, 2; 3, 4]; | |
52 ## c = condeig (a) | |
53 ## @result{} [1.0150; 1.0150] | |
54 ## @end group | |
55 ## @end example | |
56 ## @end deftypefn | |
57 | |
58 ## Author: Arno Onken <whyly@gmx.net> | |
59 ## Description: Condition numbers for eigenvalues | |
60 | |
61 function [v, lambda, c] = condeig (a) | |
62 | |
63 # Check arguments | |
64 if (nargin != 1 || nargout > 3) | |
65 usage ("[v, lambda, c] = condeig (a)"); | |
66 endif | |
67 | |
68 if (! isempty (a) && ! ismatrix (a)) | |
69 error ("condeig: a must be a numeric matrix"); | |
70 endif | |
71 | |
72 if (columns (a) != rows (a)) | |
73 error ("condeig: a must be a square matrix"); | |
74 endif | |
75 | |
76 # Right eigenvectors | |
77 [v, lambda] = eig (a); | |
78 | |
79 if (isempty (a)) | |
80 c = lambda; | |
81 else | |
82 # Corresponding left eigenvectors | |
83 vl = inv (v)'; | |
84 # Normalize vectors | |
85 vl = vl ./ repmat (sqrt (sum (abs (vl .^ 2))), rows (vl), 1); | |
86 | |
87 # Condition numbers | |
88 # cos (angle) = (norm (v1) * norm (v2)) / dot (v1, v2) | |
89 # Norm of the eigenvectors is 1 => norm (v1) * norm (v2) = 1 | |
90 c = abs (1 ./ dot (vl, v)'); | |
91 endif | |
92 | |
93 if (nargout == 0 || nargout == 1) | |
94 v = c; | |
95 endif | |
96 | |
97 endfunction | |
98 | |
99 %!test | |
100 %! a = [1, 2; 3, 4]; | |
101 %! c = condeig (a); | |
102 %! expected_c = [1.0150; 1.0150]; | |
103 %! assert (c, expected_c, 0.001); | |
104 | |
105 %!test | |
106 %! a = [1, 3; 5, 8]; | |
107 %! [v, lambda, c] = condeig (a); | |
108 %! [expected_v, expected_lambda] = eig (a); | |
109 %! expected_c = [1.0182; 1.0182]; | |
110 %! assert (v, expected_v, 0.001); | |
111 %! assert (lambda, expected_lambda, 0.001); | |
112 %! assert (c, expected_c, 0.001); |