matfuncs.py 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670
  1. #
  2. # Author: Travis Oliphant, March 2002
  3. #
  4. from __future__ import division, print_function, absolute_import
  5. __all__ = ['expm','cosm','sinm','tanm','coshm','sinhm',
  6. 'tanhm','logm','funm','signm','sqrtm',
  7. 'expm_frechet', 'expm_cond', 'fractional_matrix_power']
  8. from numpy import (Inf, dot, diag, product, logical_not, ravel,
  9. transpose, conjugate, absolute, amax, sign, isfinite, single)
  10. import numpy as np
  11. # Local imports
  12. from .misc import norm
  13. from .basic import solve, inv
  14. from .special_matrices import triu
  15. from .decomp_svd import svd
  16. from .decomp_schur import schur, rsf2csf
  17. from ._expm_frechet import expm_frechet, expm_cond
  18. from ._matfuncs_sqrtm import sqrtm
  19. eps = np.finfo(float).eps
  20. feps = np.finfo(single).eps
  21. _array_precision = {'i': 1, 'l': 1, 'f': 0, 'd': 1, 'F': 0, 'D': 1}
  22. ###############################################################################
  23. # Utility functions.
  24. def _asarray_square(A):
  25. """
  26. Wraps asarray with the extra requirement that the input be a square matrix.
  27. The motivation is that the matfuncs module has real functions that have
  28. been lifted to square matrix functions.
  29. Parameters
  30. ----------
  31. A : array_like
  32. A square matrix.
  33. Returns
  34. -------
  35. out : ndarray
  36. An ndarray copy or view or other representation of A.
  37. """
  38. A = np.asarray(A)
  39. if len(A.shape) != 2 or A.shape[0] != A.shape[1]:
  40. raise ValueError('expected square array_like input')
  41. return A
  42. def _maybe_real(A, B, tol=None):
  43. """
  44. Return either B or the real part of B, depending on properties of A and B.
  45. The motivation is that B has been computed as a complicated function of A,
  46. and B may be perturbed by negligible imaginary components.
  47. If A is real and B is complex with small imaginary components,
  48. then return a real copy of B. The assumption in that case would be that
  49. the imaginary components of B are numerical artifacts.
  50. Parameters
  51. ----------
  52. A : ndarray
  53. Input array whose type is to be checked as real vs. complex.
  54. B : ndarray
  55. Array to be returned, possibly without its imaginary part.
  56. tol : float
  57. Absolute tolerance.
  58. Returns
  59. -------
  60. out : real or complex array
  61. Either the input array B or only the real part of the input array B.
  62. """
  63. # Note that booleans and integers compare as real.
  64. if np.isrealobj(A) and np.iscomplexobj(B):
  65. if tol is None:
  66. tol = {0:feps*1e3, 1:eps*1e6}[_array_precision[B.dtype.char]]
  67. if np.allclose(B.imag, 0.0, atol=tol):
  68. B = B.real
  69. return B
  70. ###############################################################################
  71. # Matrix functions.
  72. def fractional_matrix_power(A, t):
  73. """
  74. Compute the fractional power of a matrix.
  75. Proceeds according to the discussion in section (6) of [1]_.
  76. Parameters
  77. ----------
  78. A : (N, N) array_like
  79. Matrix whose fractional power to evaluate.
  80. t : float
  81. Fractional power.
  82. Returns
  83. -------
  84. X : (N, N) array_like
  85. The fractional power of the matrix.
  86. References
  87. ----------
  88. .. [1] Nicholas J. Higham and Lijing lin (2011)
  89. "A Schur-Pade Algorithm for Fractional Powers of a Matrix."
  90. SIAM Journal on Matrix Analysis and Applications,
  91. 32 (3). pp. 1056-1078. ISSN 0895-4798
  92. Examples
  93. --------
  94. >>> from scipy.linalg import fractional_matrix_power
  95. >>> a = np.array([[1.0, 3.0], [1.0, 4.0]])
  96. >>> b = fractional_matrix_power(a, 0.5)
  97. >>> b
  98. array([[ 0.75592895, 1.13389342],
  99. [ 0.37796447, 1.88982237]])
  100. >>> np.dot(b, b) # Verify square root
  101. array([[ 1., 3.],
  102. [ 1., 4.]])
  103. """
  104. # This fixes some issue with imports;
  105. # this function calls onenormest which is in scipy.sparse.
  106. A = _asarray_square(A)
  107. import scipy.linalg._matfuncs_inv_ssq
  108. return scipy.linalg._matfuncs_inv_ssq._fractional_matrix_power(A, t)
  109. def logm(A, disp=True):
  110. """
  111. Compute matrix logarithm.
  112. The matrix logarithm is the inverse of
  113. expm: expm(logm(`A`)) == `A`
  114. Parameters
  115. ----------
  116. A : (N, N) array_like
  117. Matrix whose logarithm to evaluate
  118. disp : bool, optional
  119. Print warning if error in the result is estimated large
  120. instead of returning estimated error. (Default: True)
  121. Returns
  122. -------
  123. logm : (N, N) ndarray
  124. Matrix logarithm of `A`
  125. errest : float
  126. (if disp == False)
  127. 1-norm of the estimated error, ||err||_1 / ||A||_1
  128. References
  129. ----------
  130. .. [1] Awad H. Al-Mohy and Nicholas J. Higham (2012)
  131. "Improved Inverse Scaling and Squaring Algorithms
  132. for the Matrix Logarithm."
  133. SIAM Journal on Scientific Computing, 34 (4). C152-C169.
  134. ISSN 1095-7197
  135. .. [2] Nicholas J. Higham (2008)
  136. "Functions of Matrices: Theory and Computation"
  137. ISBN 978-0-898716-46-7
  138. .. [3] Nicholas J. Higham and Lijing lin (2011)
  139. "A Schur-Pade Algorithm for Fractional Powers of a Matrix."
  140. SIAM Journal on Matrix Analysis and Applications,
  141. 32 (3). pp. 1056-1078. ISSN 0895-4798
  142. Examples
  143. --------
  144. >>> from scipy.linalg import logm, expm
  145. >>> a = np.array([[1.0, 3.0], [1.0, 4.0]])
  146. >>> b = logm(a)
  147. >>> b
  148. array([[-1.02571087, 2.05142174],
  149. [ 0.68380725, 1.02571087]])
  150. >>> expm(b) # Verify expm(logm(a)) returns a
  151. array([[ 1., 3.],
  152. [ 1., 4.]])
  153. """
  154. A = _asarray_square(A)
  155. # Avoid circular import ... this is OK, right?
  156. import scipy.linalg._matfuncs_inv_ssq
  157. F = scipy.linalg._matfuncs_inv_ssq._logm(A)
  158. F = _maybe_real(A, F)
  159. errtol = 1000*eps
  160. #TODO use a better error approximation
  161. errest = norm(expm(F)-A,1) / norm(A,1)
  162. if disp:
  163. if not isfinite(errest) or errest >= errtol:
  164. print("logm result may be inaccurate, approximate err =", errest)
  165. return F
  166. else:
  167. return F, errest
  168. def expm(A):
  169. """
  170. Compute the matrix exponential using Pade approximation.
  171. Parameters
  172. ----------
  173. A : (N, N) array_like or sparse matrix
  174. Matrix to be exponentiated.
  175. Returns
  176. -------
  177. expm : (N, N) ndarray
  178. Matrix exponential of `A`.
  179. References
  180. ----------
  181. .. [1] Awad H. Al-Mohy and Nicholas J. Higham (2009)
  182. "A New Scaling and Squaring Algorithm for the Matrix Exponential."
  183. SIAM Journal on Matrix Analysis and Applications.
  184. 31 (3). pp. 970-989. ISSN 1095-7162
  185. Examples
  186. --------
  187. >>> from scipy.linalg import expm, sinm, cosm
  188. Matrix version of the formula exp(0) = 1:
  189. >>> expm(np.zeros((2,2)))
  190. array([[ 1., 0.],
  191. [ 0., 1.]])
  192. Euler's identity (exp(i*theta) = cos(theta) + i*sin(theta))
  193. applied to a matrix:
  194. >>> a = np.array([[1.0, 2.0], [-1.0, 3.0]])
  195. >>> expm(1j*a)
  196. array([[ 0.42645930+1.89217551j, -2.13721484-0.97811252j],
  197. [ 1.06860742+0.48905626j, -1.71075555+0.91406299j]])
  198. >>> cosm(a) + 1j*sinm(a)
  199. array([[ 0.42645930+1.89217551j, -2.13721484-0.97811252j],
  200. [ 1.06860742+0.48905626j, -1.71075555+0.91406299j]])
  201. """
  202. # Input checking and conversion is provided by sparse.linalg.expm().
  203. import scipy.sparse.linalg
  204. return scipy.sparse.linalg.expm(A)
  205. def cosm(A):
  206. """
  207. Compute the matrix cosine.
  208. This routine uses expm to compute the matrix exponentials.
  209. Parameters
  210. ----------
  211. A : (N, N) array_like
  212. Input array
  213. Returns
  214. -------
  215. cosm : (N, N) ndarray
  216. Matrix cosine of A
  217. Examples
  218. --------
  219. >>> from scipy.linalg import expm, sinm, cosm
  220. Euler's identity (exp(i*theta) = cos(theta) + i*sin(theta))
  221. applied to a matrix:
  222. >>> a = np.array([[1.0, 2.0], [-1.0, 3.0]])
  223. >>> expm(1j*a)
  224. array([[ 0.42645930+1.89217551j, -2.13721484-0.97811252j],
  225. [ 1.06860742+0.48905626j, -1.71075555+0.91406299j]])
  226. >>> cosm(a) + 1j*sinm(a)
  227. array([[ 0.42645930+1.89217551j, -2.13721484-0.97811252j],
  228. [ 1.06860742+0.48905626j, -1.71075555+0.91406299j]])
  229. """
  230. A = _asarray_square(A)
  231. if np.iscomplexobj(A):
  232. return 0.5*(expm(1j*A) + expm(-1j*A))
  233. else:
  234. return expm(1j*A).real
  235. def sinm(A):
  236. """
  237. Compute the matrix sine.
  238. This routine uses expm to compute the matrix exponentials.
  239. Parameters
  240. ----------
  241. A : (N, N) array_like
  242. Input array.
  243. Returns
  244. -------
  245. sinm : (N, N) ndarray
  246. Matrix sine of `A`
  247. Examples
  248. --------
  249. >>> from scipy.linalg import expm, sinm, cosm
  250. Euler's identity (exp(i*theta) = cos(theta) + i*sin(theta))
  251. applied to a matrix:
  252. >>> a = np.array([[1.0, 2.0], [-1.0, 3.0]])
  253. >>> expm(1j*a)
  254. array([[ 0.42645930+1.89217551j, -2.13721484-0.97811252j],
  255. [ 1.06860742+0.48905626j, -1.71075555+0.91406299j]])
  256. >>> cosm(a) + 1j*sinm(a)
  257. array([[ 0.42645930+1.89217551j, -2.13721484-0.97811252j],
  258. [ 1.06860742+0.48905626j, -1.71075555+0.91406299j]])
  259. """
  260. A = _asarray_square(A)
  261. if np.iscomplexobj(A):
  262. return -0.5j*(expm(1j*A) - expm(-1j*A))
  263. else:
  264. return expm(1j*A).imag
  265. def tanm(A):
  266. """
  267. Compute the matrix tangent.
  268. This routine uses expm to compute the matrix exponentials.
  269. Parameters
  270. ----------
  271. A : (N, N) array_like
  272. Input array.
  273. Returns
  274. -------
  275. tanm : (N, N) ndarray
  276. Matrix tangent of `A`
  277. Examples
  278. --------
  279. >>> from scipy.linalg import tanm, sinm, cosm
  280. >>> a = np.array([[1.0, 3.0], [1.0, 4.0]])
  281. >>> t = tanm(a)
  282. >>> t
  283. array([[ -2.00876993, -8.41880636],
  284. [ -2.80626879, -10.42757629]])
  285. Verify tanm(a) = sinm(a).dot(inv(cosm(a)))
  286. >>> s = sinm(a)
  287. >>> c = cosm(a)
  288. >>> s.dot(np.linalg.inv(c))
  289. array([[ -2.00876993, -8.41880636],
  290. [ -2.80626879, -10.42757629]])
  291. """
  292. A = _asarray_square(A)
  293. return _maybe_real(A, solve(cosm(A), sinm(A)))
  294. def coshm(A):
  295. """
  296. Compute the hyperbolic matrix cosine.
  297. This routine uses expm to compute the matrix exponentials.
  298. Parameters
  299. ----------
  300. A : (N, N) array_like
  301. Input array.
  302. Returns
  303. -------
  304. coshm : (N, N) ndarray
  305. Hyperbolic matrix cosine of `A`
  306. Examples
  307. --------
  308. >>> from scipy.linalg import tanhm, sinhm, coshm
  309. >>> a = np.array([[1.0, 3.0], [1.0, 4.0]])
  310. >>> c = coshm(a)
  311. >>> c
  312. array([[ 11.24592233, 38.76236492],
  313. [ 12.92078831, 50.00828725]])
  314. Verify tanhm(a) = sinhm(a).dot(inv(coshm(a)))
  315. >>> t = tanhm(a)
  316. >>> s = sinhm(a)
  317. >>> t - s.dot(np.linalg.inv(c))
  318. array([[ 2.72004641e-15, 4.55191440e-15],
  319. [ 0.00000000e+00, -5.55111512e-16]])
  320. """
  321. A = _asarray_square(A)
  322. return _maybe_real(A, 0.5 * (expm(A) + expm(-A)))
  323. def sinhm(A):
  324. """
  325. Compute the hyperbolic matrix sine.
  326. This routine uses expm to compute the matrix exponentials.
  327. Parameters
  328. ----------
  329. A : (N, N) array_like
  330. Input array.
  331. Returns
  332. -------
  333. sinhm : (N, N) ndarray
  334. Hyperbolic matrix sine of `A`
  335. Examples
  336. --------
  337. >>> from scipy.linalg import tanhm, sinhm, coshm
  338. >>> a = np.array([[1.0, 3.0], [1.0, 4.0]])
  339. >>> s = sinhm(a)
  340. >>> s
  341. array([[ 10.57300653, 39.28826594],
  342. [ 13.09608865, 49.86127247]])
  343. Verify tanhm(a) = sinhm(a).dot(inv(coshm(a)))
  344. >>> t = tanhm(a)
  345. >>> c = coshm(a)
  346. >>> t - s.dot(np.linalg.inv(c))
  347. array([[ 2.72004641e-15, 4.55191440e-15],
  348. [ 0.00000000e+00, -5.55111512e-16]])
  349. """
  350. A = _asarray_square(A)
  351. return _maybe_real(A, 0.5 * (expm(A) - expm(-A)))
  352. def tanhm(A):
  353. """
  354. Compute the hyperbolic matrix tangent.
  355. This routine uses expm to compute the matrix exponentials.
  356. Parameters
  357. ----------
  358. A : (N, N) array_like
  359. Input array
  360. Returns
  361. -------
  362. tanhm : (N, N) ndarray
  363. Hyperbolic matrix tangent of `A`
  364. Examples
  365. --------
  366. >>> from scipy.linalg import tanhm, sinhm, coshm
  367. >>> a = np.array([[1.0, 3.0], [1.0, 4.0]])
  368. >>> t = tanhm(a)
  369. >>> t
  370. array([[ 0.3428582 , 0.51987926],
  371. [ 0.17329309, 0.86273746]])
  372. Verify tanhm(a) = sinhm(a).dot(inv(coshm(a)))
  373. >>> s = sinhm(a)
  374. >>> c = coshm(a)
  375. >>> t - s.dot(np.linalg.inv(c))
  376. array([[ 2.72004641e-15, 4.55191440e-15],
  377. [ 0.00000000e+00, -5.55111512e-16]])
  378. """
  379. A = _asarray_square(A)
  380. return _maybe_real(A, solve(coshm(A), sinhm(A)))
  381. def funm(A, func, disp=True):
  382. """
  383. Evaluate a matrix function specified by a callable.
  384. Returns the value of matrix-valued function ``f`` at `A`. The
  385. function ``f`` is an extension of the scalar-valued function `func`
  386. to matrices.
  387. Parameters
  388. ----------
  389. A : (N, N) array_like
  390. Matrix at which to evaluate the function
  391. func : callable
  392. Callable object that evaluates a scalar function f.
  393. Must be vectorized (eg. using vectorize).
  394. disp : bool, optional
  395. Print warning if error in the result is estimated large
  396. instead of returning estimated error. (Default: True)
  397. Returns
  398. -------
  399. funm : (N, N) ndarray
  400. Value of the matrix function specified by func evaluated at `A`
  401. errest : float
  402. (if disp == False)
  403. 1-norm of the estimated error, ||err||_1 / ||A||_1
  404. Examples
  405. --------
  406. >>> from scipy.linalg import funm
  407. >>> a = np.array([[1.0, 3.0], [1.0, 4.0]])
  408. >>> funm(a, lambda x: x*x)
  409. array([[ 4., 15.],
  410. [ 5., 19.]])
  411. >>> a.dot(a)
  412. array([[ 4., 15.],
  413. [ 5., 19.]])
  414. Notes
  415. -----
  416. This function implements the general algorithm based on Schur decomposition
  417. (Algorithm 9.1.1. in [1]_).
  418. If the input matrix is known to be diagonalizable, then relying on the
  419. eigendecomposition is likely to be faster. For example, if your matrix is
  420. Hermitian, you can do
  421. >>> from scipy.linalg import eigh
  422. >>> def funm_herm(a, func, check_finite=False):
  423. ... w, v = eigh(a, check_finite=check_finite)
  424. ... ## if you further know that your matrix is positive semidefinite,
  425. ... ## you can optionally guard against precision errors by doing
  426. ... # w = np.maximum(w, 0)
  427. ... w = func(w)
  428. ... return (v * w).dot(v.conj().T)
  429. References
  430. ----------
  431. .. [1] Gene H. Golub, Charles F. van Loan, Matrix Computations 4th ed.
  432. """
  433. A = _asarray_square(A)
  434. # Perform Shur decomposition (lapack ?gees)
  435. T, Z = schur(A)
  436. T, Z = rsf2csf(T,Z)
  437. n,n = T.shape
  438. F = diag(func(diag(T))) # apply function to diagonal elements
  439. F = F.astype(T.dtype.char) # e.g. when F is real but T is complex
  440. minden = abs(T[0,0])
  441. # implement Algorithm 11.1.1 from Golub and Van Loan
  442. # "matrix Computations."
  443. for p in range(1,n):
  444. for i in range(1,n-p+1):
  445. j = i + p
  446. s = T[i-1,j-1] * (F[j-1,j-1] - F[i-1,i-1])
  447. ksl = slice(i,j-1)
  448. val = dot(T[i-1,ksl],F[ksl,j-1]) - dot(F[i-1,ksl],T[ksl,j-1])
  449. s = s + val
  450. den = T[j-1,j-1] - T[i-1,i-1]
  451. if den != 0.0:
  452. s = s / den
  453. F[i-1,j-1] = s
  454. minden = min(minden,abs(den))
  455. F = dot(dot(Z, F), transpose(conjugate(Z)))
  456. F = _maybe_real(A, F)
  457. tol = {0:feps, 1:eps}[_array_precision[F.dtype.char]]
  458. if minden == 0.0:
  459. minden = tol
  460. err = min(1, max(tol,(tol/minden)*norm(triu(T,1),1)))
  461. if product(ravel(logical_not(isfinite(F))),axis=0):
  462. err = Inf
  463. if disp:
  464. if err > 1000*tol:
  465. print("funm result may be inaccurate, approximate err =", err)
  466. return F
  467. else:
  468. return F, err
  469. def signm(A, disp=True):
  470. """
  471. Matrix sign function.
  472. Extension of the scalar sign(x) to matrices.
  473. Parameters
  474. ----------
  475. A : (N, N) array_like
  476. Matrix at which to evaluate the sign function
  477. disp : bool, optional
  478. Print warning if error in the result is estimated large
  479. instead of returning estimated error. (Default: True)
  480. Returns
  481. -------
  482. signm : (N, N) ndarray
  483. Value of the sign function at `A`
  484. errest : float
  485. (if disp == False)
  486. 1-norm of the estimated error, ||err||_1 / ||A||_1
  487. Examples
  488. --------
  489. >>> from scipy.linalg import signm, eigvals
  490. >>> a = [[1,2,3], [1,2,1], [1,1,1]]
  491. >>> eigvals(a)
  492. array([ 4.12488542+0.j, -0.76155718+0.j, 0.63667176+0.j])
  493. >>> eigvals(signm(a))
  494. array([-1.+0.j, 1.+0.j, 1.+0.j])
  495. """
  496. A = _asarray_square(A)
  497. def rounded_sign(x):
  498. rx = np.real(x)
  499. if rx.dtype.char == 'f':
  500. c = 1e3*feps*amax(x)
  501. else:
  502. c = 1e3*eps*amax(x)
  503. return sign((absolute(rx) > c) * rx)
  504. result, errest = funm(A, rounded_sign, disp=0)
  505. errtol = {0:1e3*feps, 1:1e3*eps}[_array_precision[result.dtype.char]]
  506. if errest < errtol:
  507. return result
  508. # Handle signm of defective matrices:
  509. # See "E.D.Denman and J.Leyva-Ramos, Appl.Math.Comp.,
  510. # 8:237-250,1981" for how to improve the following (currently a
  511. # rather naive) iteration process:
  512. # a = result # sometimes iteration converges faster but where??
  513. # Shifting to avoid zero eigenvalues. How to ensure that shifting does
  514. # not change the spectrum too much?
  515. vals = svd(A, compute_uv=0)
  516. max_sv = np.amax(vals)
  517. # min_nonzero_sv = vals[(vals>max_sv*errtol).tolist().count(1)-1]
  518. # c = 0.5/min_nonzero_sv
  519. c = 0.5/max_sv
  520. S0 = A + c*np.identity(A.shape[0])
  521. prev_errest = errest
  522. for i in range(100):
  523. iS0 = inv(S0)
  524. S0 = 0.5*(S0 + iS0)
  525. Pp = 0.5*(dot(S0,S0)+S0)
  526. errest = norm(dot(Pp,Pp)-Pp,1)
  527. if errest < errtol or prev_errest == errest:
  528. break
  529. prev_errest = errest
  530. if disp:
  531. if not isfinite(errest) or errest >= errtol:
  532. print("signm result may be inaccurate, approximate err =", errest)
  533. return S0
  534. else:
  535. return S0, errest