Mercurial > octave
annotate scripts/linear-algebra/logm.m @ 28789:28de41192f3c
Eliminate unneeded verification of nargin, nargout in m-files.
* FIRfilter.m, FIRfilter_aggregation.m, get.m, polynomial.m,
polynomial_superiorto.m, polynomial2.m, makeUniqueStrings.m, base64decode.m,
base64encode.m, cd.m, lin2mu.m, record.m, sound.m, soundsc.m, accumarray.m,
accumdim.m, bitcmp.m, bitset.m, cart2pol.m, celldisp.m, circshift.m,
cplxpair.m, cumtrapz.m, flip.m, idivide.m, interpft.m, logspace.m, pol2cart.m,
polyarea.m, postpad.m, prepad.m, rat.m, rot90.m, rotdim.m, shift.m, shiftdim.m,
sortrows.m, trapz.m, dsearch.m, dsearchn.m, getappdata.m, getpixelposition.m,
guidata.m, guihandles.m, isappdata.m, listfonts.m, uigetdir.m,
waitforbuttonpress.m, __makeinfo__.m, doc.m, get_first_help_sentence.m,
autumn.m, bone.m, brighten.m, cmpermute.m, cmunique.m, colorcube.m, contrast.m,
cool.m, copper.m, cubehelix.m, flag.m, gray.m, gray2ind.m, hot.m, hsv.m,
im2double.m, im2frame.m, imformats.m, jet.m, lines.m, ocean.m, pink.m, prism.m,
rainbow.m, rgbplot.m, spinmap.m, spring.m, summer.m, viridis.m, white.m,
winter.m, beep.m, importdata.m, is_valid_file_id.m, javachk.m, javaclasspath.m,
findstr.m, genvarname.m, strmatch.m, bandwidth.m, commutation_matrix.m, cond.m,
cross.m, isdefinite.m, ishermitian.m, issymmetric.m, krylov.m, linsolve.m,
logm.m, lscov.m, null.m, ordeig.m, orth.m, rank.m, rref.m, vecnorm.m,
bunzip2.m, citation.m, computer.m, copyfile.m, dir.m, dos.m, fileattrib.m,
gunzip.m, inputParser.m, inputname.m, ismac.m, ispc.m, isunix.m, license.m,
list_primes.m, methods.m, mkdir.m, movefile.m, nargchk.m, news.m,
orderfields.m, recycle.m, tar.m, unix.m, unpack.m, untar.m, unzip.m, ver.m,
version.m, what.m, zip.m, decic.m, fminbnd.m, fminunc.m, fsolve.m, fzero.m,
glpk.m, humps.m, lsqnonneg.m, optimget.m, pqpnonneg.m, sqp.m, pathdef.m,
camlookat.m, hidden.m, specular.m, plotmatrix.m, smooth3.m, sombrero.m,
stemleaf.m, __gnuplot_drawnow__.m, __opengl_info__.m, ancestor.m, cla.m,
close.m, closereq.m, copyobj.m, gca.m, gcf.m, ginput.m, graphics_toolkit.m,
groot.m, hgload.m, hgsave.m, isgraphics.m, ishold.m, linkaxes.m, meshgrid.m,
newplot.m, refresh.m, refreshdata.m, rotate.m, saveas.m, struct2hdl.m, conv.m,
mkpp.m, mpoles.m, padecoef.m, pchip.m, polyder.m, polyfit.m, polygcd.m,
polyint.m, polyout.m, polyval.m, ppder.m, ppint.m, getpref.m, ispref.m,
rmpref.m, profexport.m, profshow.m, powerset.m, arch_fit.m, arma_rnd.m,
blackman.m, detrend.m, diffpara.m, fftconv.m, fftfilt.m, filter2.m, freqz.m,
freqz_plot.m, hamming.m, hanning.m, sinetone.m, sinewave.m, spectral_adf.m,
spectral_xdf.m, stft.m, unwrap.m, gplot.m, ichol.m, ilu.m, spdiags.m, sprand.m,
sprandn.m, spstats.m, svds.m, treelayout.m, treeplot.m, betainc.m,
betaincinv.m, ellipke.m, gammainc.m, gammaincinv.m, legendre.m, pow2.m,
hankel.m, pascal.m, rosser.m, toeplitz.m, bounds.m, corr.m, cov.m, histc.m,
kendall.m, kurtosis.m, mad.m, mode.m, moment.m, prctile.m, quantile.m, range.m,
ranks.m, run_count.m, skewness.m, spearman.m, std.m, var.m, zscore.m,
dec2base.m, dec2bin.m, dec2hex.m, index.m, mat2str.m, native2unicode.m,
ostrsplit.m, strjoin.m, strjust.m, strtok.m, substr.m, unicode2native.m,
untabify.m, __debug_octave__.m, demo.m, example.m, fail.m, oruntests.m,
dump_demos.m, speed.m, test.m, date.m, datenum.m, datestr.m, datevec.m,
is_leap_year.m, now.m, weekday.m:
Eliminate unneeded verification of nargin, nargout in m-files now that
the interpreter checks these values.
author | Rik <rik@octave.org> |
---|---|
date | Thu, 24 Sep 2020 14:44:58 -0700 |
parents | 9f9ac219896d |
children | 90fea9cc9caa |
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 ## |
27919
1891570abac8
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
27918
diff
changeset
|
3 ## Copyright (C) 2008-2020 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
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/>. |
4334 | 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 |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
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. |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
14 ## |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
15 ## 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
|
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. |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
19 ## |
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
20 ## 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
|
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 ######################################################################## |
4334 | 25 |
26 ## -*- texinfo -*- | |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20231
diff
changeset
|
27 ## @deftypefn {} {@var{s} =} logm (@var{A}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20231
diff
changeset
|
28 ## @deftypefnx {} {@var{s} =} logm (@var{A}, @var{opt_iters}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20231
diff
changeset
|
29 ## @deftypefnx {} {[@var{s}, @var{iters}] =} logm (@dots{}) |
20160
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
30 ## Compute the matrix logarithm of the square matrix @var{A}. |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
31 ## |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
32 ## The 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
|
33 ## |
10827
228cd18455a6
logm.m: Improve documentation string. Add GPL header. Add additional test block.
Rik <octave@nomad.inbox5.com>
parents:
10825
diff
changeset
|
34 ## @example |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
35 ## 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
|
36 ## @end example |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
37 ## |
20160
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
38 ## The optional input @var{opt_iters} is the maximum number of square roots |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
39 ## to compute and defaults to 100. |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
40 ## |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
41 ## The optional output @var{iters} is the number of square roots actually |
03b9d17a2d95
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
42 ## computed. |
12584
7ef7e20057fa
Improve documentation strings in Linear Algebra chapter.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
43 ## @seealso{expm, sqrtm} |
4334 | 44 ## @end deftypefn |
45 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
46 ## 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
|
47 ## (SIAM, 2008.) |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
48 ## |
4334 | 49 |
16768 | 50 ## Author: N. J. Higham |
51 ## Author: Richard T. Guy <guyrt7@wfu.edu> | |
27985
9f9ac219896d
maint: Remove remaining "Author:" instances from code base.
Rik <rik@octave.org>
parents:
27978
diff
changeset
|
52 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
53 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
|
54 |
28789
28de41192f3c
Eliminate unneeded verification of nargin, nargout in m-files.
Rik <rik@octave.org>
parents:
27985
diff
changeset
|
55 if (nargin == 0) |
6046 | 56 print_usage (); |
4334 | 57 endif |
58 | |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
59 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
|
60 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
|
61 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
62 |
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
|
63 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
|
64 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
|
65 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
|
66 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
|
67 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
|
68 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
|
69 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
|
70 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
71 [u, s] = schur (A); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
72 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
73 if (isreal (A)) |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
74 [u, s] = rsf2csf (u, s); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
75 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
76 |
13097
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
77 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
|
78 if (any (eigv < 0)) |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
79 warning ("Octave:logm:non-principal", |
11458 | 80 "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
|
81 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
82 |
13097
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
83 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
|
84 |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
85 k = 0; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
86 ## 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
|
87 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
|
88 p = 0; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
89 m = 7; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
90 while (k < opt_iters) |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
91 tau = norm (s - eye (size (s)),1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
92 if (tau <= theta (7)) |
20231
83792dd9bcc1
Use in-place operators in m-files where possible.
Rik <rik@octave.org>
parents:
20160
diff
changeset
|
93 p += 1; |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
94 j(1) = find (tau <= theta, 1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
95 j(2) = find (tau / 2 <= theta, 1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
96 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
|
97 m = j(1); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
98 break |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
99 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
100 endif |
20231
83792dd9bcc1
Use in-place operators in m-files where possible.
Rik <rik@octave.org>
parents:
20160
diff
changeset
|
101 k += 1; |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
102 s = sqrtm (s); |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
103 endwhile |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
104 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
105 if (k >= opt_iters) |
11458 | 106 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
|
107 endif |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
108 |
20231
83792dd9bcc1
Use in-place operators in m-files where possible.
Rik <rik@octave.org>
parents:
20160
diff
changeset
|
109 s -= eye (size (s)); |
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
|
110 |
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
|
111 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
|
112 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
|
113 endif |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
114 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
115 s = 2^k * u * s * u'; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
116 |
13097
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
117 ## Remove small complex values (O(eps)) which may have entered calculation |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14491
diff
changeset
|
118 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
|
119 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
|
120 endif |
52e4aa30d5b2
logm.m: Return real matrix when all eigenvalues are real (Bug #32121).
Rik <octave@nomad.inbox5.com>
parents:
12584
diff
changeset
|
121 |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
122 if (nargout == 2) |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
123 iters = k; |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
124 endif |
4334 | 125 |
126 endfunction | |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
127 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
128 ################## ANCILLARY FUNCTIONS ################################ |
11033
d9c8916bb9dd
Untabify a few remaining .m scripts.
Rik <octave@nomad.inbox5.com>
parents:
10827
diff
changeset
|
129 ###### 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
|
130 ###### Reference: |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
131 ###### 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
|
132 ###### (SIAM, 2008.). |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
133 ####################################################################### |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
134 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
135 ##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
|
136 ## 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
|
137 ## 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
|
138 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
139 function s = logm_pade_pf (A, m) |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
140 |
11458 | 141 [nodes, wts] = gauss_legendre (m); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
142 ## Convert from [-1,1] to [0,1]. |
11458 | 143 nodes = (nodes+1)/2; |
20231
83792dd9bcc1
Use in-place operators in m-files where possible.
Rik <rik@octave.org>
parents:
20160
diff
changeset
|
144 wts /= 2; |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
145 |
11471
994e2a93a8e2
Use uppercase 'A' to refer to matrix inputs in m-files.
Rik <octave@nomad.inbox5.com>
parents:
11465
diff
changeset
|
146 n = length (A); |
11458 | 147 s = zeros (n); |
148 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
|
149 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
|
150 endfor |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
151 |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
152 endfunction |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
153 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
154 ###################################################################### |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
155 ## 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
|
156 ## [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
|
157 ## for N-point Gauss-Legendre quadrature. |
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 ## Reference: |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
160 ## 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
|
161 ## 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
|
162 |
11458 | 163 function [x, w] = gauss_legendre (n) |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
164 |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
165 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
|
166 v = i./sqrt ((2*i).^2-1); |
11458 | 167 [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
|
168 x = diag (D); |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
169 w = 2*(V(1,:)'.^2); |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
170 |
10825
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
171 endfunction |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
172 |
cace99cb01ab
rewrite logm (M. Caliari, R.T. Guy)
Jaroslav Hajek <highegg@gmail.com>
parents:
7017
diff
changeset
|
173 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
174 %!assert (norm (logm ([1 -1;0 1]) - [0 -1; 0 0]) < 1e-5) |
19670
e5facc6eec13
Silence warning messages in %!test code.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
175 %!test |
e5facc6eec13
Silence warning messages in %!test code.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
176 %! warning ("off", "Octave:logm:non-principal", "local"); |
e5facc6eec13
Silence warning messages in %!test code.
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
177 %! assert (norm (expm (logm ([-1 2 ; 4 -1])) - [-1 2 ; 4 -1]) < 1e-5); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
178 %!assert (logm ([1 -1 -1;0 1 -1; 0 0 1]), [0 -1 -1.5; 0 0 -1; 0 0 0], 1e-5) |
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
|
179 %!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
|
180 %!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
|
181 %!assert (full (logm (10*eye (3))), logm (full (10*eye (3))), 8*eps) |
14481
d2bffa78730e
Fix logm for complex matrix with real eigenvalues (bug #34893).
Marco Caliari <marco.caliari@univr.it>
parents:
14327
diff
changeset
|
182 %!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
|
183 |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
184 ## Test input validation |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
185 %!error logm () |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
186 %!error logm (1, 2, 3) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
187 %!error <logm: A must be a square matrix> logm ([1 0;0 1; 2 2]) |