File: logm.m

package info (click to toggle)
octave 9.4.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 144,300 kB
  • sloc: cpp: 332,784; ansic: 77,239; fortran: 20,963; objc: 9,396; sh: 8,213; yacc: 4,925; lex: 4,389; perl: 1,544; java: 1,366; awk: 1,259; makefile: 648; xml: 189
file content (211 lines) | stat: -rw-r--r-- 6,310 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
########################################################################
##
## Copyright (C) 2008-2024 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## This file is part of Octave.
##
## Octave is free software: you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation, either version 3 of the License, or
## (at your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
## GNU General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <https://www.gnu.org/licenses/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn  {} {@var{s} =} logm (@var{A})
## @deftypefnx {} {@var{s} =} logm (@var{A}, @var{opt_iters})
## @deftypefnx {} {[@var{s}, @var{iters}] =} logm (@dots{})
## Compute the matrix logarithm of the square matrix @var{A}.
##
## The implementation utilizes a Pad@'e approximant and the identity
##
## @example
## logm (@var{A}) = 2^k * logm (@var{A}^(1 / 2^k))
## @end example
##
## The optional input @var{opt_iters} is the maximum number of square roots
## to compute and defaults to 100.
##
## The optional output @var{iters} is the number of square roots actually
## computed.
## @seealso{expm, sqrtm}
## @end deftypefn

## Reference: N. J. Higham, Functions of Matrices: Theory and Computation
##            (SIAM, 2008.)
##

## Author: N. J. Higham
## Author: Richard T. Guy <guyrt7@wfu.edu>

function [s, iters] = logm (A, opt_iters = 100)

  if (nargin == 0)
    print_usage ();
  endif

  if (! issquare (A))
    error ("logm: A must be a square matrix");
  endif

  if (isscalar (A))
    s = log (A);
    return;
  elseif (isdiag (A))
    s = diag (log (diag (A)));
    return;
  endif

  [u, s] = schur (A);

  if (isreal (A))
    [u, s] = rsf2csf (u, s);
  endif

  eigv = diag (s);
  n = rows (A);
  tol = n * eps (max (abs (eigv)));
  real_neg_eigv = (real (eigv) < -tol) & (imag (eigv) <= tol);
  if (any (real_neg_eigv))
    warning ("Octave:logm:non-principal",
             "logm: principal matrix logarithm is not defined for matrices with negative eigenvalues; computing non-principal logarithm");
  endif

  real_eig = ! any (real_neg_eigv);

  if (max (abs (triu (s,1))(:)) < tol)
    ## Will run for Hermitian matrices as Schur decomposition is diagonal.
    ## This way is faster and more accurate but only works on a diagonal matrix.
    logeigv = log (eigv);
    logeigv(isinf (logeigv)) = -log (realmax ());
    s = u * diag (logeigv) * u';
    iters = 0;
  else
    k = 0;
    ## Algorithm 11.9 in "Function of matrices", by N. Higham
    theta = [0, 0, 1.61e-2, 5.38e-2, 1.13e-1, 1.86e-1, 2.6429608311114350e-1];
    p = 0;
    m = 7;
    while (k < opt_iters)
      tau = norm (s - eye (n), 1);
      if (tau <= theta (7))
        p += 1;
        j(1) = find (tau <= theta, 1);
        j(2) = find (tau / 2 <= theta, 1);
        if (j(1) - j(2) <= 1 || p == 2)
          m = j(1);
          break;
        endif
      endif
      k += 1;
      s = sqrtm (s);
    endwhile

    if (k >= opt_iters)
      warning ("logm: maximum number of square roots exceeded; results may still be accurate");
    endif

    s -= eye (n);

    if (m > 1)
      s = logm_pade_pf (s, m);
    endif

    s = 2^k * u * s * u';

    if (nargout == 2)
      iters = k;
    endif
  endif
  ## Remove small complex values (O(eps)) which may have entered calculation
  if (real_eig && isreal (A))
    s = real (s);
  endif

endfunction

################## ANCILLARY FUNCTIONS ################################
######  Taken from the mfttoolbox (GPL 3) by D. Higham.
######  Reference:
######      D. Higham, Functions of Matrices: Theory and Computation
######      (SIAM, 2008.).
#######################################################################

## LOGM_PADE_PF   Evaluate Pade approximant to matrix log by partial fractions.
##   Y = LOGM_PADE_PF(A,M) evaluates the [M/M] Pade approximation to
##   LOG(EYE(SIZE(A))+A) using a partial fraction expansion.

function s = logm_pade_pf (A, m)

  [nodes, wts] = gauss_legendre (m);
  ## Convert from [-1,1] to [0,1].
  nodes = (nodes+1)/2;
  wts /= 2;

  n = length (A);
  s = zeros (n);
  for j = 1:m
    s += wts(j)*(A/(eye (n) + nodes(j)*A));
  endfor

endfunction

######################################################################
## GAUSS_LEGENDRE  Nodes and weights for Gauss-Legendre quadrature.
##   [X,W] = GAUSS_LEGENDRE(N) computes the nodes X and weights W
##   for N-point Gauss-Legendre quadrature.

## Reference:
## G. H. Golub and J. H. Welsch, Calculation of Gauss quadrature
## rules, Math. Comp., 23(106):221-230, 1969.

function [x, w] = gauss_legendre (n)

  i = 1:n-1;
  v = i./sqrt ((2*i).^2-1);
  [V, D] = eig (diag (v, -1) + diag (v, 1));
  x = diag (D);
  w = 2*(V(1,:)'.^2);

endfunction


%!assert (norm (logm ([1 -1;0 1]) - [0 -1; 0 0]) < 1e-5)
%!test
%! warning ("off", "Octave:logm:non-principal", "local");
%! assert (norm (expm (logm ([-1 2 ; 4 -1])) - [-1 2 ; 4 -1]) < 1e-5);
%!assert (logm ([1 -1 -1;0 1 -1; 0 0 1]), [0 -1 -1.5; 0 0 -1; 0 0 0], 1e-5)
%!assert (logm (10), log (10))
%!assert (full (logm (eye (3))), logm (full (eye (3))))
%!assert (full (logm (10*eye (3))), logm (full (10*eye (3))), 8*eps)
%!assert (logm (expm ([0 1i; -1i 0])), [0 1i; -1i 0], 10 * eps)
%!test <*60738>
%! A = [0.2510, 1.2808, -1.2252; ...
%!      0.2015, 1.0766, 0.5630; ...
%!      -1.9769, -1.0922, -0.5831];
%! if (__have_feature__ ("LLVM_LIBCXX"))
%!   ## The math libraries in libc++ seem to require larger tolerances
%!   tol = 65*eps;
%! else
%!   tol = 40*eps;
%! endif
%! warning ("off", "Octave:logm:non-principal", "local");
%! assert (expm (logm (A)), A, tol);
%!assert (expm (logm (eye (3))), eye (3))
%!assert (expm (logm (zeros (3))), zeros (3))

## Test input validation
%!error <Invalid call> logm ()
%!error <logm: A must be a square matrix> logm ([1 0;0 1; 2 2])