_expm_frechet.py 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411
  1. """Frechet derivative of the matrix exponential."""
  2. from __future__ import division, print_function, absolute_import
  3. import numpy as np
  4. import scipy.linalg
  5. __all__ = ['expm_frechet', 'expm_cond']
  6. def expm_frechet(A, E, method=None, compute_expm=True, check_finite=True):
  7. """
  8. Frechet derivative of the matrix exponential of A in the direction E.
  9. Parameters
  10. ----------
  11. A : (N, N) array_like
  12. Matrix of which to take the matrix exponential.
  13. E : (N, N) array_like
  14. Matrix direction in which to take the Frechet derivative.
  15. method : str, optional
  16. Choice of algorithm. Should be one of
  17. - `SPS` (default)
  18. - `blockEnlarge`
  19. compute_expm : bool, optional
  20. Whether to compute also `expm_A` in addition to `expm_frechet_AE`.
  21. Default is True.
  22. check_finite : bool, optional
  23. Whether to check that the input matrix contains only finite numbers.
  24. Disabling may give a performance gain, but may result in problems
  25. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  26. Returns
  27. -------
  28. expm_A : ndarray
  29. Matrix exponential of A.
  30. expm_frechet_AE : ndarray
  31. Frechet derivative of the matrix exponential of A in the direction E.
  32. For ``compute_expm = False``, only `expm_frechet_AE` is returned.
  33. See also
  34. --------
  35. expm : Compute the exponential of a matrix.
  36. Notes
  37. -----
  38. This section describes the available implementations that can be selected
  39. by the `method` parameter. The default method is *SPS*.
  40. Method *blockEnlarge* is a naive algorithm.
  41. Method *SPS* is Scaling-Pade-Squaring [1]_.
  42. It is a sophisticated implementation which should take
  43. only about 3/8 as much time as the naive implementation.
  44. The asymptotics are the same.
  45. .. versionadded:: 0.13.0
  46. References
  47. ----------
  48. .. [1] Awad H. Al-Mohy and Nicholas J. Higham (2009)
  49. Computing the Frechet Derivative of the Matrix Exponential,
  50. with an application to Condition Number Estimation.
  51. SIAM Journal On Matrix Analysis and Applications.,
  52. 30 (4). pp. 1639-1657. ISSN 1095-7162
  53. Examples
  54. --------
  55. >>> import scipy.linalg
  56. >>> A = np.random.randn(3, 3)
  57. >>> E = np.random.randn(3, 3)
  58. >>> expm_A, expm_frechet_AE = scipy.linalg.expm_frechet(A, E)
  59. >>> expm_A.shape, expm_frechet_AE.shape
  60. ((3, 3), (3, 3))
  61. >>> import scipy.linalg
  62. >>> A = np.random.randn(3, 3)
  63. >>> E = np.random.randn(3, 3)
  64. >>> expm_A, expm_frechet_AE = scipy.linalg.expm_frechet(A, E)
  65. >>> M = np.zeros((6, 6))
  66. >>> M[:3, :3] = A; M[:3, 3:] = E; M[3:, 3:] = A
  67. >>> expm_M = scipy.linalg.expm(M)
  68. >>> np.allclose(expm_A, expm_M[:3, :3])
  69. True
  70. >>> np.allclose(expm_frechet_AE, expm_M[:3, 3:])
  71. True
  72. """
  73. if check_finite:
  74. A = np.asarray_chkfinite(A)
  75. E = np.asarray_chkfinite(E)
  76. else:
  77. A = np.asarray(A)
  78. E = np.asarray(E)
  79. if A.ndim != 2 or A.shape[0] != A.shape[1]:
  80. raise ValueError('expected A to be a square matrix')
  81. if E.ndim != 2 or E.shape[0] != E.shape[1]:
  82. raise ValueError('expected E to be a square matrix')
  83. if A.shape != E.shape:
  84. raise ValueError('expected A and E to be the same shape')
  85. if method is None:
  86. method = 'SPS'
  87. if method == 'SPS':
  88. expm_A, expm_frechet_AE = expm_frechet_algo_64(A, E)
  89. elif method == 'blockEnlarge':
  90. expm_A, expm_frechet_AE = expm_frechet_block_enlarge(A, E)
  91. else:
  92. raise ValueError('Unknown implementation %s' % method)
  93. if compute_expm:
  94. return expm_A, expm_frechet_AE
  95. else:
  96. return expm_frechet_AE
  97. def expm_frechet_block_enlarge(A, E):
  98. """
  99. This is a helper function, mostly for testing and profiling.
  100. Return expm(A), frechet(A, E)
  101. """
  102. n = A.shape[0]
  103. M = np.vstack([
  104. np.hstack([A, E]),
  105. np.hstack([np.zeros_like(A), A])])
  106. expm_M = scipy.linalg.expm(M)
  107. return expm_M[:n, :n], expm_M[:n, n:]
  108. """
  109. Maximal values ell_m of ||2**-s A|| such that the backward error bound
  110. does not exceed 2**-53.
  111. """
  112. ell_table_61 = (
  113. None,
  114. # 1
  115. 2.11e-8,
  116. 3.56e-4,
  117. 1.08e-2,
  118. 6.49e-2,
  119. 2.00e-1,
  120. 4.37e-1,
  121. 7.83e-1,
  122. 1.23e0,
  123. 1.78e0,
  124. 2.42e0,
  125. # 11
  126. 3.13e0,
  127. 3.90e0,
  128. 4.74e0,
  129. 5.63e0,
  130. 6.56e0,
  131. 7.52e0,
  132. 8.53e0,
  133. 9.56e0,
  134. 1.06e1,
  135. 1.17e1,
  136. )
  137. # The b vectors and U and V are copypasted
  138. # from scipy.sparse.linalg.matfuncs.py.
  139. # M, Lu, Lv follow (6.11), (6.12), (6.13), (3.3)
  140. def _diff_pade3(A, E, ident):
  141. b = (120., 60., 12., 1.)
  142. A2 = A.dot(A)
  143. M2 = np.dot(A, E) + np.dot(E, A)
  144. U = A.dot(b[3]*A2 + b[1]*ident)
  145. V = b[2]*A2 + b[0]*ident
  146. Lu = A.dot(b[3]*M2) + E.dot(b[3]*A2 + b[1]*ident)
  147. Lv = b[2]*M2
  148. return U, V, Lu, Lv
  149. def _diff_pade5(A, E, ident):
  150. b = (30240., 15120., 3360., 420., 30., 1.)
  151. A2 = A.dot(A)
  152. M2 = np.dot(A, E) + np.dot(E, A)
  153. A4 = np.dot(A2, A2)
  154. M4 = np.dot(A2, M2) + np.dot(M2, A2)
  155. U = A.dot(b[5]*A4 + b[3]*A2 + b[1]*ident)
  156. V = b[4]*A4 + b[2]*A2 + b[0]*ident
  157. Lu = (A.dot(b[5]*M4 + b[3]*M2) +
  158. E.dot(b[5]*A4 + b[3]*A2 + b[1]*ident))
  159. Lv = b[4]*M4 + b[2]*M2
  160. return U, V, Lu, Lv
  161. def _diff_pade7(A, E, ident):
  162. b = (17297280., 8648640., 1995840., 277200., 25200., 1512., 56., 1.)
  163. A2 = A.dot(A)
  164. M2 = np.dot(A, E) + np.dot(E, A)
  165. A4 = np.dot(A2, A2)
  166. M4 = np.dot(A2, M2) + np.dot(M2, A2)
  167. A6 = np.dot(A2, A4)
  168. M6 = np.dot(A4, M2) + np.dot(M4, A2)
  169. U = A.dot(b[7]*A6 + b[5]*A4 + b[3]*A2 + b[1]*ident)
  170. V = b[6]*A6 + b[4]*A4 + b[2]*A2 + b[0]*ident
  171. Lu = (A.dot(b[7]*M6 + b[5]*M4 + b[3]*M2) +
  172. E.dot(b[7]*A6 + b[5]*A4 + b[3]*A2 + b[1]*ident))
  173. Lv = b[6]*M6 + b[4]*M4 + b[2]*M2
  174. return U, V, Lu, Lv
  175. def _diff_pade9(A, E, ident):
  176. b = (17643225600., 8821612800., 2075673600., 302702400., 30270240.,
  177. 2162160., 110880., 3960., 90., 1.)
  178. A2 = A.dot(A)
  179. M2 = np.dot(A, E) + np.dot(E, A)
  180. A4 = np.dot(A2, A2)
  181. M4 = np.dot(A2, M2) + np.dot(M2, A2)
  182. A6 = np.dot(A2, A4)
  183. M6 = np.dot(A4, M2) + np.dot(M4, A2)
  184. A8 = np.dot(A4, A4)
  185. M8 = np.dot(A4, M4) + np.dot(M4, A4)
  186. U = A.dot(b[9]*A8 + b[7]*A6 + b[5]*A4 + b[3]*A2 + b[1]*ident)
  187. V = b[8]*A8 + b[6]*A6 + b[4]*A4 + b[2]*A2 + b[0]*ident
  188. Lu = (A.dot(b[9]*M8 + b[7]*M6 + b[5]*M4 + b[3]*M2) +
  189. E.dot(b[9]*A8 + b[7]*A6 + b[5]*A4 + b[3]*A2 + b[1]*ident))
  190. Lv = b[8]*M8 + b[6]*M6 + b[4]*M4 + b[2]*M2
  191. return U, V, Lu, Lv
  192. def expm_frechet_algo_64(A, E):
  193. n = A.shape[0]
  194. s = None
  195. ident = np.identity(n)
  196. A_norm_1 = scipy.linalg.norm(A, 1)
  197. m_pade_pairs = (
  198. (3, _diff_pade3),
  199. (5, _diff_pade5),
  200. (7, _diff_pade7),
  201. (9, _diff_pade9))
  202. for m, pade in m_pade_pairs:
  203. if A_norm_1 <= ell_table_61[m]:
  204. U, V, Lu, Lv = pade(A, E, ident)
  205. s = 0
  206. break
  207. if s is None:
  208. # scaling
  209. s = max(0, int(np.ceil(np.log2(A_norm_1 / ell_table_61[13]))))
  210. A = A * 2.0**-s
  211. E = E * 2.0**-s
  212. # pade order 13
  213. A2 = np.dot(A, A)
  214. M2 = np.dot(A, E) + np.dot(E, A)
  215. A4 = np.dot(A2, A2)
  216. M4 = np.dot(A2, M2) + np.dot(M2, A2)
  217. A6 = np.dot(A2, A4)
  218. M6 = np.dot(A4, M2) + np.dot(M4, A2)
  219. b = (64764752532480000., 32382376266240000., 7771770303897600.,
  220. 1187353796428800., 129060195264000., 10559470521600.,
  221. 670442572800., 33522128640., 1323241920., 40840800., 960960.,
  222. 16380., 182., 1.)
  223. W1 = b[13]*A6 + b[11]*A4 + b[9]*A2
  224. W2 = b[7]*A6 + b[5]*A4 + b[3]*A2 + b[1]*ident
  225. Z1 = b[12]*A6 + b[10]*A4 + b[8]*A2
  226. Z2 = b[6]*A6 + b[4]*A4 + b[2]*A2 + b[0]*ident
  227. W = np.dot(A6, W1) + W2
  228. U = np.dot(A, W)
  229. V = np.dot(A6, Z1) + Z2
  230. Lw1 = b[13]*M6 + b[11]*M4 + b[9]*M2
  231. Lw2 = b[7]*M6 + b[5]*M4 + b[3]*M2
  232. Lz1 = b[12]*M6 + b[10]*M4 + b[8]*M2
  233. Lz2 = b[6]*M6 + b[4]*M4 + b[2]*M2
  234. Lw = np.dot(A6, Lw1) + np.dot(M6, W1) + Lw2
  235. Lu = np.dot(A, Lw) + np.dot(E, W)
  236. Lv = np.dot(A6, Lz1) + np.dot(M6, Z1) + Lz2
  237. # factor once and solve twice
  238. lu_piv = scipy.linalg.lu_factor(-U + V)
  239. R = scipy.linalg.lu_solve(lu_piv, U + V)
  240. L = scipy.linalg.lu_solve(lu_piv, Lu + Lv + np.dot((Lu - Lv), R))
  241. # squaring
  242. for k in range(s):
  243. L = np.dot(R, L) + np.dot(L, R)
  244. R = np.dot(R, R)
  245. return R, L
  246. def vec(M):
  247. """
  248. Stack columns of M to construct a single vector.
  249. This is somewhat standard notation in linear algebra.
  250. Parameters
  251. ----------
  252. M : 2d array_like
  253. Input matrix
  254. Returns
  255. -------
  256. v : 1d ndarray
  257. Output vector
  258. """
  259. return M.T.ravel()
  260. def expm_frechet_kronform(A, method=None, check_finite=True):
  261. """
  262. Construct the Kronecker form of the Frechet derivative of expm.
  263. Parameters
  264. ----------
  265. A : array_like with shape (N, N)
  266. Matrix to be expm'd.
  267. method : str, optional
  268. Extra keyword to be passed to expm_frechet.
  269. check_finite : bool, optional
  270. Whether to check that the input matrix contains only finite numbers.
  271. Disabling may give a performance gain, but may result in problems
  272. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  273. Returns
  274. -------
  275. K : 2d ndarray with shape (N*N, N*N)
  276. Kronecker form of the Frechet derivative of the matrix exponential.
  277. Notes
  278. -----
  279. This function is used to help compute the condition number
  280. of the matrix exponential.
  281. See also
  282. --------
  283. expm : Compute a matrix exponential.
  284. expm_frechet : Compute the Frechet derivative of the matrix exponential.
  285. expm_cond : Compute the relative condition number of the matrix exponential
  286. in the Frobenius norm.
  287. """
  288. if check_finite:
  289. A = np.asarray_chkfinite(A)
  290. else:
  291. A = np.asarray(A)
  292. if len(A.shape) != 2 or A.shape[0] != A.shape[1]:
  293. raise ValueError('expected a square matrix')
  294. n = A.shape[0]
  295. ident = np.identity(n)
  296. cols = []
  297. for i in range(n):
  298. for j in range(n):
  299. E = np.outer(ident[i], ident[j])
  300. F = expm_frechet(A, E,
  301. method=method, compute_expm=False, check_finite=False)
  302. cols.append(vec(F))
  303. return np.vstack(cols).T
  304. def expm_cond(A, check_finite=True):
  305. """
  306. Relative condition number of the matrix exponential in the Frobenius norm.
  307. Parameters
  308. ----------
  309. A : 2d array_like
  310. Square input matrix with shape (N, N).
  311. check_finite : bool, optional
  312. Whether to check that the input matrix contains only finite numbers.
  313. Disabling may give a performance gain, but may result in problems
  314. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  315. Returns
  316. -------
  317. kappa : float
  318. The relative condition number of the matrix exponential
  319. in the Frobenius norm
  320. Notes
  321. -----
  322. A faster estimate for the condition number in the 1-norm
  323. has been published but is not yet implemented in scipy.
  324. .. versionadded:: 0.14.0
  325. See also
  326. --------
  327. expm : Compute the exponential of a matrix.
  328. expm_frechet : Compute the Frechet derivative of the matrix exponential.
  329. Examples
  330. --------
  331. >>> from scipy.linalg import expm_cond
  332. >>> A = np.array([[-0.3, 0.2, 0.6], [0.6, 0.3, -0.1], [-0.7, 1.2, 0.9]])
  333. >>> k = expm_cond(A)
  334. >>> k
  335. 1.7787805864469866
  336. """
  337. if check_finite:
  338. A = np.asarray_chkfinite(A)
  339. else:
  340. A = np.asarray(A)
  341. if len(A.shape) != 2 or A.shape[0] != A.shape[1]:
  342. raise ValueError('expected a square matrix')
  343. X = scipy.linalg.expm(A)
  344. K = expm_frechet_kronform(A, check_finite=False)
  345. # The following norm choices are deliberate.
  346. # The norms of A and X are Frobenius norms,
  347. # and the norm of K is the induced 2-norm.
  348. A_norm = scipy.linalg.norm(A, 'fro')
  349. X_norm = scipy.linalg.norm(X, 'fro')
  350. K_norm = scipy.linalg.norm(K, 2)
  351. kappa = (K_norm * A_norm) / X_norm
  352. return kappa