matfuncs.py 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868
  1. """
  2. Sparse matrix functions
  3. """
  4. #
  5. # Authors: Travis Oliphant, March 2002
  6. # Anthony Scopatz, August 2012 (Sparse Updates)
  7. # Jake Vanderplas, August 2012 (Sparse Updates)
  8. #
  9. from __future__ import division, print_function, absolute_import
  10. __all__ = ['expm', 'inv']
  11. import math
  12. import numpy as np
  13. import scipy.special
  14. from scipy.linalg.basic import solve, solve_triangular
  15. from scipy.sparse.base import isspmatrix
  16. from scipy.sparse.construct import eye as speye
  17. from scipy.sparse.linalg import spsolve
  18. import scipy.sparse
  19. import scipy.sparse.linalg
  20. from scipy.sparse.linalg.interface import LinearOperator
  21. UPPER_TRIANGULAR = 'upper_triangular'
  22. def inv(A):
  23. """
  24. Compute the inverse of a sparse matrix
  25. Parameters
  26. ----------
  27. A : (M,M) ndarray or sparse matrix
  28. square matrix to be inverted
  29. Returns
  30. -------
  31. Ainv : (M,M) ndarray or sparse matrix
  32. inverse of `A`
  33. Notes
  34. -----
  35. This computes the sparse inverse of `A`. If the inverse of `A` is expected
  36. to be non-sparse, it will likely be faster to convert `A` to dense and use
  37. scipy.linalg.inv.
  38. Examples
  39. --------
  40. >>> from scipy.sparse import csc_matrix
  41. >>> from scipy.sparse.linalg import inv
  42. >>> A = csc_matrix([[1., 0.], [1., 2.]])
  43. >>> Ainv = inv(A)
  44. >>> Ainv
  45. <2x2 sparse matrix of type '<class 'numpy.float64'>'
  46. with 3 stored elements in Compressed Sparse Column format>
  47. >>> A.dot(Ainv)
  48. <2x2 sparse matrix of type '<class 'numpy.float64'>'
  49. with 2 stored elements in Compressed Sparse Column format>
  50. >>> A.dot(Ainv).todense()
  51. matrix([[ 1., 0.],
  52. [ 0., 1.]])
  53. .. versionadded:: 0.12.0
  54. """
  55. #check input
  56. if not scipy.sparse.isspmatrix(A):
  57. raise TypeError('Input must be a sparse matrix')
  58. I = speye(A.shape[0], A.shape[1], dtype=A.dtype, format=A.format)
  59. Ainv = spsolve(A, I)
  60. return Ainv
  61. def _onenorm_matrix_power_nnm(A, p):
  62. """
  63. Compute the 1-norm of a non-negative integer power of a non-negative matrix.
  64. Parameters
  65. ----------
  66. A : a square ndarray or matrix or sparse matrix
  67. Input matrix with non-negative entries.
  68. p : non-negative integer
  69. The power to which the matrix is to be raised.
  70. Returns
  71. -------
  72. out : float
  73. The 1-norm of the matrix power p of A.
  74. """
  75. # check input
  76. if int(p) != p or p < 0:
  77. raise ValueError('expected non-negative integer p')
  78. p = int(p)
  79. if len(A.shape) != 2 or A.shape[0] != A.shape[1]:
  80. raise ValueError('expected A to be like a square matrix')
  81. # Explicitly make a column vector so that this works when A is a
  82. # numpy matrix (in addition to ndarray and sparse matrix).
  83. v = np.ones((A.shape[0], 1), dtype=float)
  84. M = A.T
  85. for i in range(p):
  86. v = M.dot(v)
  87. return np.max(v)
  88. def _onenorm(A):
  89. # A compatibility function which should eventually disappear.
  90. # This is copypasted from expm_action.
  91. if scipy.sparse.isspmatrix(A):
  92. return max(abs(A).sum(axis=0).flat)
  93. else:
  94. return np.linalg.norm(A, 1)
  95. def _ident_like(A):
  96. # A compatibility function which should eventually disappear.
  97. # This is copypasted from expm_action.
  98. if scipy.sparse.isspmatrix(A):
  99. return scipy.sparse.construct.eye(A.shape[0], A.shape[1],
  100. dtype=A.dtype, format=A.format)
  101. else:
  102. return np.eye(A.shape[0], A.shape[1], dtype=A.dtype)
  103. def _is_upper_triangular(A):
  104. # This function could possibly be of wider interest.
  105. if isspmatrix(A):
  106. lower_part = scipy.sparse.tril(A, -1)
  107. # Check structural upper triangularity,
  108. # then coincidental upper triangularity if needed.
  109. return lower_part.nnz == 0 or lower_part.count_nonzero() == 0
  110. else:
  111. return not np.tril(A, -1).any()
  112. def _smart_matrix_product(A, B, alpha=None, structure=None):
  113. """
  114. A matrix product that knows about sparse and structured matrices.
  115. Parameters
  116. ----------
  117. A : 2d ndarray
  118. First matrix.
  119. B : 2d ndarray
  120. Second matrix.
  121. alpha : float
  122. The matrix product will be scaled by this constant.
  123. structure : str, optional
  124. A string describing the structure of both matrices `A` and `B`.
  125. Only `upper_triangular` is currently supported.
  126. Returns
  127. -------
  128. M : 2d ndarray
  129. Matrix product of A and B.
  130. """
  131. if len(A.shape) != 2:
  132. raise ValueError('expected A to be a rectangular matrix')
  133. if len(B.shape) != 2:
  134. raise ValueError('expected B to be a rectangular matrix')
  135. f = None
  136. if structure == UPPER_TRIANGULAR:
  137. if not isspmatrix(A) and not isspmatrix(B):
  138. f, = scipy.linalg.get_blas_funcs(('trmm',), (A, B))
  139. if f is not None:
  140. if alpha is None:
  141. alpha = 1.
  142. out = f(alpha, A, B)
  143. else:
  144. if alpha is None:
  145. out = A.dot(B)
  146. else:
  147. out = alpha * A.dot(B)
  148. return out
  149. class MatrixPowerOperator(LinearOperator):
  150. def __init__(self, A, p, structure=None):
  151. if A.ndim != 2 or A.shape[0] != A.shape[1]:
  152. raise ValueError('expected A to be like a square matrix')
  153. if p < 0:
  154. raise ValueError('expected p to be a non-negative integer')
  155. self._A = A
  156. self._p = p
  157. self._structure = structure
  158. self.dtype = A.dtype
  159. self.ndim = A.ndim
  160. self.shape = A.shape
  161. def _matvec(self, x):
  162. for i in range(self._p):
  163. x = self._A.dot(x)
  164. return x
  165. def _rmatvec(self, x):
  166. A_T = self._A.T
  167. x = x.ravel()
  168. for i in range(self._p):
  169. x = A_T.dot(x)
  170. return x
  171. def _matmat(self, X):
  172. for i in range(self._p):
  173. X = _smart_matrix_product(self._A, X, structure=self._structure)
  174. return X
  175. @property
  176. def T(self):
  177. return MatrixPowerOperator(self._A.T, self._p)
  178. class ProductOperator(LinearOperator):
  179. """
  180. For now, this is limited to products of multiple square matrices.
  181. """
  182. def __init__(self, *args, **kwargs):
  183. self._structure = kwargs.get('structure', None)
  184. for A in args:
  185. if len(A.shape) != 2 or A.shape[0] != A.shape[1]:
  186. raise ValueError(
  187. 'For now, the ProductOperator implementation is '
  188. 'limited to the product of multiple square matrices.')
  189. if args:
  190. n = args[0].shape[0]
  191. for A in args:
  192. for d in A.shape:
  193. if d != n:
  194. raise ValueError(
  195. 'The square matrices of the ProductOperator '
  196. 'must all have the same shape.')
  197. self.shape = (n, n)
  198. self.ndim = len(self.shape)
  199. self.dtype = np.find_common_type([x.dtype for x in args], [])
  200. self._operator_sequence = args
  201. def _matvec(self, x):
  202. for A in reversed(self._operator_sequence):
  203. x = A.dot(x)
  204. return x
  205. def _rmatvec(self, x):
  206. x = x.ravel()
  207. for A in self._operator_sequence:
  208. x = A.T.dot(x)
  209. return x
  210. def _matmat(self, X):
  211. for A in reversed(self._operator_sequence):
  212. X = _smart_matrix_product(A, X, structure=self._structure)
  213. return X
  214. @property
  215. def T(self):
  216. T_args = [A.T for A in reversed(self._operator_sequence)]
  217. return ProductOperator(*T_args)
  218. def _onenormest_matrix_power(A, p,
  219. t=2, itmax=5, compute_v=False, compute_w=False, structure=None):
  220. """
  221. Efficiently estimate the 1-norm of A^p.
  222. Parameters
  223. ----------
  224. A : ndarray
  225. Matrix whose 1-norm of a power is to be computed.
  226. p : int
  227. Non-negative integer power.
  228. t : int, optional
  229. A positive parameter controlling the tradeoff between
  230. accuracy versus time and memory usage.
  231. Larger values take longer and use more memory
  232. but give more accurate output.
  233. itmax : int, optional
  234. Use at most this many iterations.
  235. compute_v : bool, optional
  236. Request a norm-maximizing linear operator input vector if True.
  237. compute_w : bool, optional
  238. Request a norm-maximizing linear operator output vector if True.
  239. Returns
  240. -------
  241. est : float
  242. An underestimate of the 1-norm of the sparse matrix.
  243. v : ndarray, optional
  244. The vector such that ||Av||_1 == est*||v||_1.
  245. It can be thought of as an input to the linear operator
  246. that gives an output with particularly large norm.
  247. w : ndarray, optional
  248. The vector Av which has relatively large 1-norm.
  249. It can be thought of as an output of the linear operator
  250. that is relatively large in norm compared to the input.
  251. """
  252. return scipy.sparse.linalg.onenormest(
  253. MatrixPowerOperator(A, p, structure=structure))
  254. def _onenormest_product(operator_seq,
  255. t=2, itmax=5, compute_v=False, compute_w=False, structure=None):
  256. """
  257. Efficiently estimate the 1-norm of the matrix product of the args.
  258. Parameters
  259. ----------
  260. operator_seq : linear operator sequence
  261. Matrices whose 1-norm of product is to be computed.
  262. t : int, optional
  263. A positive parameter controlling the tradeoff between
  264. accuracy versus time and memory usage.
  265. Larger values take longer and use more memory
  266. but give more accurate output.
  267. itmax : int, optional
  268. Use at most this many iterations.
  269. compute_v : bool, optional
  270. Request a norm-maximizing linear operator input vector if True.
  271. compute_w : bool, optional
  272. Request a norm-maximizing linear operator output vector if True.
  273. structure : str, optional
  274. A string describing the structure of all operators.
  275. Only `upper_triangular` is currently supported.
  276. Returns
  277. -------
  278. est : float
  279. An underestimate of the 1-norm of the sparse matrix.
  280. v : ndarray, optional
  281. The vector such that ||Av||_1 == est*||v||_1.
  282. It can be thought of as an input to the linear operator
  283. that gives an output with particularly large norm.
  284. w : ndarray, optional
  285. The vector Av which has relatively large 1-norm.
  286. It can be thought of as an output of the linear operator
  287. that is relatively large in norm compared to the input.
  288. """
  289. return scipy.sparse.linalg.onenormest(
  290. ProductOperator(*operator_seq, structure=structure))
  291. class _ExpmPadeHelper(object):
  292. """
  293. Help lazily evaluate a matrix exponential.
  294. The idea is to not do more work than we need for high expm precision,
  295. so we lazily compute matrix powers and store or precompute
  296. other properties of the matrix.
  297. """
  298. def __init__(self, A, structure=None, use_exact_onenorm=False):
  299. """
  300. Initialize the object.
  301. Parameters
  302. ----------
  303. A : a dense or sparse square numpy matrix or ndarray
  304. The matrix to be exponentiated.
  305. structure : str, optional
  306. A string describing the structure of matrix `A`.
  307. Only `upper_triangular` is currently supported.
  308. use_exact_onenorm : bool, optional
  309. If True then only the exact one-norm of matrix powers and products
  310. will be used. Otherwise, the one-norm of powers and products
  311. may initially be estimated.
  312. """
  313. self.A = A
  314. self._A2 = None
  315. self._A4 = None
  316. self._A6 = None
  317. self._A8 = None
  318. self._A10 = None
  319. self._d4_exact = None
  320. self._d6_exact = None
  321. self._d8_exact = None
  322. self._d10_exact = None
  323. self._d4_approx = None
  324. self._d6_approx = None
  325. self._d8_approx = None
  326. self._d10_approx = None
  327. self.ident = _ident_like(A)
  328. self.structure = structure
  329. self.use_exact_onenorm = use_exact_onenorm
  330. @property
  331. def A2(self):
  332. if self._A2 is None:
  333. self._A2 = _smart_matrix_product(
  334. self.A, self.A, structure=self.structure)
  335. return self._A2
  336. @property
  337. def A4(self):
  338. if self._A4 is None:
  339. self._A4 = _smart_matrix_product(
  340. self.A2, self.A2, structure=self.structure)
  341. return self._A4
  342. @property
  343. def A6(self):
  344. if self._A6 is None:
  345. self._A6 = _smart_matrix_product(
  346. self.A4, self.A2, structure=self.structure)
  347. return self._A6
  348. @property
  349. def A8(self):
  350. if self._A8 is None:
  351. self._A8 = _smart_matrix_product(
  352. self.A6, self.A2, structure=self.structure)
  353. return self._A8
  354. @property
  355. def A10(self):
  356. if self._A10 is None:
  357. self._A10 = _smart_matrix_product(
  358. self.A4, self.A6, structure=self.structure)
  359. return self._A10
  360. @property
  361. def d4_tight(self):
  362. if self._d4_exact is None:
  363. self._d4_exact = _onenorm(self.A4)**(1/4.)
  364. return self._d4_exact
  365. @property
  366. def d6_tight(self):
  367. if self._d6_exact is None:
  368. self._d6_exact = _onenorm(self.A6)**(1/6.)
  369. return self._d6_exact
  370. @property
  371. def d8_tight(self):
  372. if self._d8_exact is None:
  373. self._d8_exact = _onenorm(self.A8)**(1/8.)
  374. return self._d8_exact
  375. @property
  376. def d10_tight(self):
  377. if self._d10_exact is None:
  378. self._d10_exact = _onenorm(self.A10)**(1/10.)
  379. return self._d10_exact
  380. @property
  381. def d4_loose(self):
  382. if self.use_exact_onenorm:
  383. return self.d4_tight
  384. if self._d4_exact is not None:
  385. return self._d4_exact
  386. else:
  387. if self._d4_approx is None:
  388. self._d4_approx = _onenormest_matrix_power(self.A2, 2,
  389. structure=self.structure)**(1/4.)
  390. return self._d4_approx
  391. @property
  392. def d6_loose(self):
  393. if self.use_exact_onenorm:
  394. return self.d6_tight
  395. if self._d6_exact is not None:
  396. return self._d6_exact
  397. else:
  398. if self._d6_approx is None:
  399. self._d6_approx = _onenormest_matrix_power(self.A2, 3,
  400. structure=self.structure)**(1/6.)
  401. return self._d6_approx
  402. @property
  403. def d8_loose(self):
  404. if self.use_exact_onenorm:
  405. return self.d8_tight
  406. if self._d8_exact is not None:
  407. return self._d8_exact
  408. else:
  409. if self._d8_approx is None:
  410. self._d8_approx = _onenormest_matrix_power(self.A4, 2,
  411. structure=self.structure)**(1/8.)
  412. return self._d8_approx
  413. @property
  414. def d10_loose(self):
  415. if self.use_exact_onenorm:
  416. return self.d10_tight
  417. if self._d10_exact is not None:
  418. return self._d10_exact
  419. else:
  420. if self._d10_approx is None:
  421. self._d10_approx = _onenormest_product((self.A4, self.A6),
  422. structure=self.structure)**(1/10.)
  423. return self._d10_approx
  424. def pade3(self):
  425. b = (120., 60., 12., 1.)
  426. U = _smart_matrix_product(self.A,
  427. b[3]*self.A2 + b[1]*self.ident,
  428. structure=self.structure)
  429. V = b[2]*self.A2 + b[0]*self.ident
  430. return U, V
  431. def pade5(self):
  432. b = (30240., 15120., 3360., 420., 30., 1.)
  433. U = _smart_matrix_product(self.A,
  434. b[5]*self.A4 + b[3]*self.A2 + b[1]*self.ident,
  435. structure=self.structure)
  436. V = b[4]*self.A4 + b[2]*self.A2 + b[0]*self.ident
  437. return U, V
  438. def pade7(self):
  439. b = (17297280., 8648640., 1995840., 277200., 25200., 1512., 56., 1.)
  440. U = _smart_matrix_product(self.A,
  441. b[7]*self.A6 + b[5]*self.A4 + b[3]*self.A2 + b[1]*self.ident,
  442. structure=self.structure)
  443. V = b[6]*self.A6 + b[4]*self.A4 + b[2]*self.A2 + b[0]*self.ident
  444. return U, V
  445. def pade9(self):
  446. b = (17643225600., 8821612800., 2075673600., 302702400., 30270240.,
  447. 2162160., 110880., 3960., 90., 1.)
  448. U = _smart_matrix_product(self.A,
  449. (b[9]*self.A8 + b[7]*self.A6 + b[5]*self.A4 +
  450. b[3]*self.A2 + b[1]*self.ident),
  451. structure=self.structure)
  452. V = (b[8]*self.A8 + b[6]*self.A6 + b[4]*self.A4 +
  453. b[2]*self.A2 + b[0]*self.ident)
  454. return U, V
  455. def pade13_scaled(self, s):
  456. b = (64764752532480000., 32382376266240000., 7771770303897600.,
  457. 1187353796428800., 129060195264000., 10559470521600.,
  458. 670442572800., 33522128640., 1323241920., 40840800., 960960.,
  459. 16380., 182., 1.)
  460. B = self.A * 2**-s
  461. B2 = self.A2 * 2**(-2*s)
  462. B4 = self.A4 * 2**(-4*s)
  463. B6 = self.A6 * 2**(-6*s)
  464. U2 = _smart_matrix_product(B6,
  465. b[13]*B6 + b[11]*B4 + b[9]*B2,
  466. structure=self.structure)
  467. U = _smart_matrix_product(B,
  468. (U2 + b[7]*B6 + b[5]*B4 +
  469. b[3]*B2 + b[1]*self.ident),
  470. structure=self.structure)
  471. V2 = _smart_matrix_product(B6,
  472. b[12]*B6 + b[10]*B4 + b[8]*B2,
  473. structure=self.structure)
  474. V = V2 + b[6]*B6 + b[4]*B4 + b[2]*B2 + b[0]*self.ident
  475. return U, V
  476. def expm(A):
  477. """
  478. Compute the matrix exponential using Pade approximation.
  479. Parameters
  480. ----------
  481. A : (M,M) array_like or sparse matrix
  482. 2D Array or Matrix (sparse or dense) to be exponentiated
  483. Returns
  484. -------
  485. expA : (M,M) ndarray
  486. Matrix exponential of `A`
  487. Notes
  488. -----
  489. This is algorithm (6.1) which is a simplification of algorithm (5.1).
  490. .. versionadded:: 0.12.0
  491. References
  492. ----------
  493. .. [1] Awad H. Al-Mohy and Nicholas J. Higham (2009)
  494. "A New Scaling and Squaring Algorithm for the Matrix Exponential."
  495. SIAM Journal on Matrix Analysis and Applications.
  496. 31 (3). pp. 970-989. ISSN 1095-7162
  497. Examples
  498. --------
  499. >>> from scipy.sparse import csc_matrix
  500. >>> from scipy.sparse.linalg import expm
  501. >>> A = csc_matrix([[1, 0, 0], [0, 2, 0], [0, 0, 3]])
  502. >>> A.todense()
  503. matrix([[1, 0, 0],
  504. [0, 2, 0],
  505. [0, 0, 3]], dtype=int64)
  506. >>> Aexp = expm(A)
  507. >>> Aexp
  508. <3x3 sparse matrix of type '<class 'numpy.float64'>'
  509. with 3 stored elements in Compressed Sparse Column format>
  510. >>> Aexp.todense()
  511. matrix([[ 2.71828183, 0. , 0. ],
  512. [ 0. , 7.3890561 , 0. ],
  513. [ 0. , 0. , 20.08553692]])
  514. """
  515. return _expm(A, use_exact_onenorm='auto')
  516. def _expm(A, use_exact_onenorm):
  517. # Core of expm, separated to allow testing exact and approximate
  518. # algorithms.
  519. # Avoid indiscriminate asarray() to allow sparse or other strange arrays.
  520. if isinstance(A, (list, tuple, np.matrix)):
  521. A = np.asarray(A)
  522. if len(A.shape) != 2 or A.shape[0] != A.shape[1]:
  523. raise ValueError('expected a square matrix')
  524. # Trivial case
  525. if A.shape == (1, 1):
  526. out = [[np.exp(A[0, 0])]]
  527. # Avoid indiscriminate casting to ndarray to
  528. # allow for sparse or other strange arrays
  529. if isspmatrix(A):
  530. return A.__class__(out)
  531. return np.array(out)
  532. # Ensure input is of float type, to avoid integer overflows etc.
  533. if ((isinstance(A, np.ndarray) or isspmatrix(A))
  534. and not np.issubdtype(A.dtype, np.inexact)):
  535. A = A.astype(float)
  536. # Detect upper triangularity.
  537. structure = UPPER_TRIANGULAR if _is_upper_triangular(A) else None
  538. if use_exact_onenorm == "auto":
  539. # Hardcode a matrix order threshold for exact vs. estimated one-norms.
  540. use_exact_onenorm = A.shape[0] < 200
  541. # Track functions of A to help compute the matrix exponential.
  542. h = _ExpmPadeHelper(
  543. A, structure=structure, use_exact_onenorm=use_exact_onenorm)
  544. # Try Pade order 3.
  545. eta_1 = max(h.d4_loose, h.d6_loose)
  546. if eta_1 < 1.495585217958292e-002 and _ell(h.A, 3) == 0:
  547. U, V = h.pade3()
  548. return _solve_P_Q(U, V, structure=structure)
  549. # Try Pade order 5.
  550. eta_2 = max(h.d4_tight, h.d6_loose)
  551. if eta_2 < 2.539398330063230e-001 and _ell(h.A, 5) == 0:
  552. U, V = h.pade5()
  553. return _solve_P_Q(U, V, structure=structure)
  554. # Try Pade orders 7 and 9.
  555. eta_3 = max(h.d6_tight, h.d8_loose)
  556. if eta_3 < 9.504178996162932e-001 and _ell(h.A, 7) == 0:
  557. U, V = h.pade7()
  558. return _solve_P_Q(U, V, structure=structure)
  559. if eta_3 < 2.097847961257068e+000 and _ell(h.A, 9) == 0:
  560. U, V = h.pade9()
  561. return _solve_P_Q(U, V, structure=structure)
  562. # Use Pade order 13.
  563. eta_4 = max(h.d8_loose, h.d10_loose)
  564. eta_5 = min(eta_3, eta_4)
  565. theta_13 = 4.25
  566. # Choose smallest s>=0 such that 2**(-s) eta_5 <= theta_13
  567. if eta_5 == 0:
  568. # Nilpotent special case
  569. s = 0
  570. else:
  571. s = max(int(np.ceil(np.log2(eta_5 / theta_13))), 0)
  572. s = s + _ell(2**-s * h.A, 13)
  573. U, V = h.pade13_scaled(s)
  574. X = _solve_P_Q(U, V, structure=structure)
  575. if structure == UPPER_TRIANGULAR:
  576. # Invoke Code Fragment 2.1.
  577. X = _fragment_2_1(X, h.A, s)
  578. else:
  579. # X = r_13(A)^(2^s) by repeated squaring.
  580. for i in range(s):
  581. X = X.dot(X)
  582. return X
  583. def _solve_P_Q(U, V, structure=None):
  584. """
  585. A helper function for expm_2009.
  586. Parameters
  587. ----------
  588. U : ndarray
  589. Pade numerator.
  590. V : ndarray
  591. Pade denominator.
  592. structure : str, optional
  593. A string describing the structure of both matrices `U` and `V`.
  594. Only `upper_triangular` is currently supported.
  595. Notes
  596. -----
  597. The `structure` argument is inspired by similar args
  598. for theano and cvxopt functions.
  599. """
  600. P = U + V
  601. Q = -U + V
  602. if isspmatrix(U):
  603. return spsolve(Q, P)
  604. elif structure is None:
  605. return solve(Q, P)
  606. elif structure == UPPER_TRIANGULAR:
  607. return solve_triangular(Q, P)
  608. else:
  609. raise ValueError('unsupported matrix structure: ' + str(structure))
  610. def _sinch(x):
  611. """
  612. Stably evaluate sinch.
  613. Notes
  614. -----
  615. The strategy of falling back to a sixth order Taylor expansion
  616. was suggested by the Spallation Neutron Source docs
  617. which was found on the internet by google search.
  618. http://www.ornl.gov/~t6p/resources/xal/javadoc/gov/sns/tools/math/ElementaryFunction.html
  619. The details of the cutoff point and the Horner-like evaluation
  620. was picked without reference to anything in particular.
  621. Note that sinch is not currently implemented in scipy.special,
  622. whereas the "engineer's" definition of sinc is implemented.
  623. The implementation of sinc involves a scaling factor of pi
  624. that distinguishes it from the "mathematician's" version of sinc.
  625. """
  626. # If x is small then use sixth order Taylor expansion.
  627. # How small is small? I am using the point where the relative error
  628. # of the approximation is less than 1e-14.
  629. # If x is large then directly evaluate sinh(x) / x.
  630. x2 = x*x
  631. if abs(x) < 0.0135:
  632. return 1 + (x2/6.)*(1 + (x2/20.)*(1 + (x2/42.)))
  633. else:
  634. return np.sinh(x) / x
  635. def _eq_10_42(lam_1, lam_2, t_12):
  636. """
  637. Equation (10.42) of Functions of Matrices: Theory and Computation.
  638. Notes
  639. -----
  640. This is a helper function for _fragment_2_1 of expm_2009.
  641. Equation (10.42) is on page 251 in the section on Schur algorithms.
  642. In particular, section 10.4.3 explains the Schur-Parlett algorithm.
  643. expm([[lam_1, t_12], [0, lam_1])
  644. =
  645. [[exp(lam_1), t_12*exp((lam_1 + lam_2)/2)*sinch((lam_1 - lam_2)/2)],
  646. [0, exp(lam_2)]
  647. """
  648. # The plain formula t_12 * (exp(lam_2) - exp(lam_2)) / (lam_2 - lam_1)
  649. # apparently suffers from cancellation, according to Higham's textbook.
  650. # A nice implementation of sinch, defined as sinh(x)/x,
  651. # will apparently work around the cancellation.
  652. a = 0.5 * (lam_1 + lam_2)
  653. b = 0.5 * (lam_1 - lam_2)
  654. return t_12 * np.exp(a) * _sinch(b)
  655. def _fragment_2_1(X, T, s):
  656. """
  657. A helper function for expm_2009.
  658. Notes
  659. -----
  660. The argument X is modified in-place, but this modification is not the same
  661. as the returned value of the function.
  662. This function also takes pains to do things in ways that are compatible
  663. with sparse matrices, for example by avoiding fancy indexing
  664. and by using methods of the matrices whenever possible instead of
  665. using functions of the numpy or scipy libraries themselves.
  666. """
  667. # Form X = r_m(2^-s T)
  668. # Replace diag(X) by exp(2^-s diag(T)).
  669. n = X.shape[0]
  670. diag_T = np.ravel(T.diagonal().copy())
  671. # Replace diag(X) by exp(2^-s diag(T)).
  672. scale = 2 ** -s
  673. exp_diag = np.exp(scale * diag_T)
  674. for k in range(n):
  675. X[k, k] = exp_diag[k]
  676. for i in range(s-1, -1, -1):
  677. X = X.dot(X)
  678. # Replace diag(X) by exp(2^-i diag(T)).
  679. scale = 2 ** -i
  680. exp_diag = np.exp(scale * diag_T)
  681. for k in range(n):
  682. X[k, k] = exp_diag[k]
  683. # Replace (first) superdiagonal of X by explicit formula
  684. # for superdiagonal of exp(2^-i T) from Eq (10.42) of
  685. # the author's 2008 textbook
  686. # Functions of Matrices: Theory and Computation.
  687. for k in range(n-1):
  688. lam_1 = scale * diag_T[k]
  689. lam_2 = scale * diag_T[k+1]
  690. t_12 = scale * T[k, k+1]
  691. value = _eq_10_42(lam_1, lam_2, t_12)
  692. X[k, k+1] = value
  693. # Return the updated X matrix.
  694. return X
  695. def _ell(A, m):
  696. """
  697. A helper function for expm_2009.
  698. Parameters
  699. ----------
  700. A : linear operator
  701. A linear operator whose norm of power we care about.
  702. m : int
  703. The power of the linear operator
  704. Returns
  705. -------
  706. value : int
  707. A value related to a bound.
  708. """
  709. if len(A.shape) != 2 or A.shape[0] != A.shape[1]:
  710. raise ValueError('expected A to be like a square matrix')
  711. p = 2*m + 1
  712. # The c_i are explained in (2.2) and (2.6) of the 2005 expm paper.
  713. # They are coefficients of terms of a generating function series expansion.
  714. choose_2p_p = scipy.special.comb(2*p, p, exact=True)
  715. abs_c_recip = float(choose_2p_p * math.factorial(2*p + 1))
  716. # This is explained after Eq. (1.2) of the 2009 expm paper.
  717. # It is the "unit roundoff" of IEEE double precision arithmetic.
  718. u = 2**-53
  719. # Compute the one-norm of matrix power p of abs(A).
  720. A_abs_onenorm = _onenorm_matrix_power_nnm(abs(A), p)
  721. # Treat zero norm as a special case.
  722. if not A_abs_onenorm:
  723. return 0
  724. alpha = A_abs_onenorm / (_onenorm(A) * abs_c_recip)
  725. log2_alpha_div_u = np.log2(alpha/u)
  726. value = int(np.ceil(log2_alpha_div_u / (2 * m)))
  727. return max(value, 0)