Mercurial > hg > octave-lyh
annotate scripts/linear-algebra/logm.m @ 14481:d2bffa78730e stable
Fix logm for complex matrix with real eigenvalues (bug #34893).
* crsf2csf, zrsf2csf: Fix off-by-one error.
* logm.m: Only truncate imaginary parts for real matrices. Add a test.
* schur.cc: Add a test for rsf2csf.x
author | Marco Caliari <marco.caliari@univr.it> |
---|---|
date | Tue, 13 Mar 2012 11:56:35 +0100 |
parents | 4d917a6a858b |
children | 5bd9e47e9277 |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
13097
diff
changeset
|
1 ## Copyright (C) 2008-2012 N.J. Higham |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
2 ## Copyright (C) 2010 Richard T. Guy <guyrt7@wfu.edu> |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
3 ## Copyright (C) 2010 Marco Caliari <marco.caliari@univr.it> |
4334 | 4 ## |
5 ## This file is part of Octave. | |
6 ## | |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
7 ## Octave is free software; you can redistribute it and/or modify it |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
8 ## under the terms of the GNU General Public License as published by |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
9 ## the Free Software Foundation; either version 3 of the License, or (at |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
10 ## your option) any later version. |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
11 ## |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
12 ## Octave is distributed in the hope that it will be useful, but |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
13 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
14 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
15 ## General Public License for more details. |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
16 ## |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
17 ## You should have received a copy of the GNU General Public License |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
18 ## along with Octave; see the file COPYING. If not, see |
7016 | 19 ## <http://www.gnu.org/licenses/>. |
4334 | 20 |
21 ## -*- texinfo -*- | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
22 ## @deftypefn {Function File} {@var{s} =} logm (@var{A}) |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
23 ## @deftypefnx {Function File} {@var{s} =} logm (@var{A}, @var{opt_iters}) |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
24 ## @deftypefnx {Function File} {[@var{s}, @var{iters}] =} logm (@dots{}) |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
25 ## Compute the matrix logarithm of the square matrix @var{A}. The |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
26 ## implementation utilizes a Pad@'e approximant and the identity |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
27 ## |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
28 ## @example |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
29 ## logm (@var{A}) = 2^k * logm (@var{A}^(1 / 2^k)) |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
30 ## @end example |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
31 ## |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
32 ## The optional argument @var{opt_iters} is the maximum number of square roots |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
33 ## to compute and defaults to 100. The optional output @var{iters} is the |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
34 ## number of square roots actually computed. |
12584
7ef7e20057fa
Improve documentation strings in Linear Algebra chapter.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
35 ## @seealso{expm, sqrtm} |
4334 | 36 ## @end deftypefn |
37 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
38 ## Reference: N. J. Higham, Functions of Matrices: Theory and Computation |
11033
d9c8916bb9dd
Untabify a few remaining .m scripts.
Rik <octave@nomad.inbox5.com>
parents:
10827
diff
changeset
|
39 ## (SIAM, 2008.) |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
40 ## |
4334 | 41 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
42 function [s, iters] = logm (A, opt_iters = 100) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
43 |
11458 | 44 if (nargin == 0 || nargin > 2) |
6046 | 45 print_usage (); |
4334 | 46 endif |
47 | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
48 if (! issquare (A)) |
11472
1740012184f9
Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents:
11471
diff
changeset
|
49 error ("logm: A must be a square matrix"); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
50 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
51 |
11476
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
52 if (isscalar (A)) |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
53 s = log (A); |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
54 return; |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
55 elseif (strfind (typeinfo (A), "diagonal matrix")) |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
56 s = diag (log (diag (A))); |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
57 return; |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
58 endif |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
59 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
60 [u, s] = schur (A); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
61 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
62 if (isreal (A)) |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
63 [u, s] = rsf2csf (u, s); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
64 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
65 |
13097
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
66 eigv = diag (s); |
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
67 if (any (eigv < 0)) |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
68 warning ("Octave:logm:non-principal", |
11458 | 69 "logm: principal matrix logarithm is not defined for matrices with negative eigenvalues; computing non-principal logarithm"); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
70 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
71 |
13097
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
72 real_eig = all (eigv >= 0); |
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
73 |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
74 k = 0; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
75 ## Algorithm 11.9 in "Function of matrices", by N. Higham |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
76 theta = [0, 0, 1.61e-2, 5.38e-2, 1.13e-1, 1.86e-1, 2.6429608311114350e-1]; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
77 p = 0; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
78 m = 7; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
79 while (k < opt_iters) |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
80 tau = norm (s - eye (size (s)),1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
81 if (tau <= theta (7)) |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
82 p = p + 1; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
83 j(1) = find (tau <= theta, 1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
84 j(2) = find (tau / 2 <= theta, 1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
85 if (j(1) - j(2) <= 1 || p == 2) |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
86 m = j(1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
87 break |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
88 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
89 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
90 k = k + 1; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
91 s = sqrtm (s); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
92 endwhile |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
93 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
94 if (k >= opt_iters) |
11458 | 95 warning ("logm: maximum number of square roots exceeded; results may still be accurate"); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
96 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
97 |
11476
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
98 s = s - eye (size (s)); |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
99 |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
100 if (m > 1) |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
101 s = logm_pade_pf (s, m); |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
102 endif |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
103 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
104 s = 2^k * u * s * u'; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
105 |
13097
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
106 ## Remove small complex values (O(eps)) which may have entered calculation |
14481
d2bffa78730e
Fix logm for complex matrix with real eigenvalues (bug #34893).
Marco Caliari <marco.caliari@univr.it>
parents:
14327
diff
changeset
|
107 if (real_eig && isreal(A)) |
13097
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
108 s = real (s); |
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
109 endif |
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
110 |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
111 if (nargout == 2) |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
112 iters = k; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
113 endif |
4334 | 114 |
115 endfunction | |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
116 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
117 ################## ANCILLARY FUNCTIONS ################################ |
11033
d9c8916bb9dd
Untabify a few remaining .m scripts.
Rik <octave@nomad.inbox5.com>
parents:
10827
diff
changeset
|
118 ###### Taken from the mfttoolbox (GPL 3) by D. Higham. |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
119 ###### Reference: |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
120 ###### D. Higham, Functions of Matrices: Theory and Computation |
11033
d9c8916bb9dd
Untabify a few remaining .m scripts.
Rik <octave@nomad.inbox5.com>
parents:
10827
diff
changeset
|
121 ###### (SIAM, 2008.). |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
122 ####################################################################### |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
123 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
124 ##LOGM_PADE_PF Evaluate Pade approximant to matrix log by partial fractions. |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
125 ## Y = LOGM_PADE_PF(A,M) evaluates the [M/M] Pade approximation to |
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
126 ## LOG(EYE(SIZE(A))+A) using a partial fraction expansion. |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
127 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
128 function s = logm_pade_pf (A, m) |
11458 | 129 [nodes, wts] = gauss_legendre (m); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
130 ## Convert from [-1,1] to [0,1]. |
11458 | 131 nodes = (nodes+1)/2; |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
132 wts = wts/2; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
133 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
134 n = length (A); |
11458 | 135 s = zeros (n); |
136 for j = 1:m | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
137 s += wts(j)*(A/(eye (n) + nodes(j)*A)); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
138 endfor |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
139 endfunction |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
140 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
141 ###################################################################### |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
142 ## GAUSS_LEGENDRE Nodes and weights for Gauss-Legendre quadrature. |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
143 ## [X,W] = GAUSS_LEGENDRE(N) computes the nodes X and weights W |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
144 ## for N-point Gauss-Legendre quadrature. |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
145 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
146 ## Reference: |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
147 ## G. H. Golub and J. H. Welsch, Calculation of Gauss quadrature |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
148 ## rules, Math. Comp., 23(106):221-230, 1969. |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
149 |
11458 | 150 function [x, w] = gauss_legendre (n) |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
151 i = 1:n-1; |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
152 v = i./sqrt ((2*i).^2-1); |
11458 | 153 [V, D] = eig (diag (v, -1) + diag (v, 1)); |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
154 x = diag (D); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
155 w = 2*(V(1,:)'.^2); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
156 endfunction |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
157 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
158 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
159 %!assert(norm(logm([1 -1;0 1]) - [0 -1; 0 0]) < 1e-5); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
160 %!assert(norm(expm(logm([-1 2 ; 4 -1])) - [-1 2 ; 4 -1]) < 1e-5); |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
161 %!assert(logm([1 -1 -1;0 1 -1; 0 0 1]), [0 -1 -1.5; 0 0 -1; 0 0 0], 1e-5); |
14481
d2bffa78730e
Fix logm for complex matrix with real eigenvalues (bug #34893).
Marco Caliari <marco.caliari@univr.it>
parents:
14327
diff
changeset
|
162 %!assert (logm (expm ([0 1i; -1i 0])), [0 1i; -1i 0], 10 * eps) |
11472
1740012184f9
Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents:
11471
diff
changeset
|
163 |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
164 %% Test input validation |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
165 %!error logm (); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
166 %!error logm (1, 2, 3); |
11472
1740012184f9
Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents:
11471
diff
changeset
|
167 %!error <logm: A must be a square matrix> logm([1 0;0 1; 2 2]); |
11476
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
168 |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
169 %!assert (logm (10), log (10)) |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
170 %!assert (full (logm (eye (3))), logm (full (eye (3)))) |
ff7e0776ba0f
logm.m: handle scalar and diagonal matrix arguments specially; call log_pade_pf only if m > 1
John W. Eaton <jwe@octave.org>
parents:
11472
diff
changeset
|
171 %!assert (full (logm (10*eye (3))), logm (full (10*eye (3))), 8*eps) |