eigen.py 40 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365
  1. from types import FunctionType
  2. from collections import Counter
  3. from mpmath import mp, workprec
  4. from mpmath.libmp.libmpf import prec_to_dps
  5. from sympy.core.sorting import default_sort_key
  6. from sympy.core.evalf import DEFAULT_MAXPREC, PrecisionExhausted
  7. from sympy.core.logic import fuzzy_and, fuzzy_or
  8. from sympy.core.numbers import Float
  9. from sympy.core.sympify import _sympify
  10. from sympy.functions.elementary.miscellaneous import sqrt
  11. from sympy.polys import roots, CRootOf, ZZ, QQ, EX
  12. from sympy.polys.matrices import DomainMatrix
  13. from sympy.polys.matrices.eigen import dom_eigenvects, dom_eigenvects_to_sympy
  14. from sympy.simplify import nsimplify, simplify as _simplify
  15. from sympy.utilities.exceptions import sympy_deprecation_warning
  16. from .common import MatrixError, NonSquareMatrixError
  17. from .determinant import _find_reasonable_pivot
  18. from .utilities import _iszero
  19. def _eigenvals_eigenvects_mpmath(M):
  20. norm2 = lambda v: mp.sqrt(sum(i**2 for i in v))
  21. v1 = None
  22. prec = max([x._prec for x in M.atoms(Float)])
  23. eps = 2**-prec
  24. while prec < DEFAULT_MAXPREC:
  25. with workprec(prec):
  26. A = mp.matrix(M.evalf(n=prec_to_dps(prec)))
  27. E, ER = mp.eig(A)
  28. v2 = norm2([i for e in E for i in (mp.re(e), mp.im(e))])
  29. if v1 is not None and mp.fabs(v1 - v2) < eps:
  30. return E, ER
  31. v1 = v2
  32. prec *= 2
  33. # we get here because the next step would have taken us
  34. # past MAXPREC or because we never took a step; in case
  35. # of the latter, we refuse to send back a solution since
  36. # it would not have been verified; we also resist taking
  37. # a small step to arrive exactly at MAXPREC since then
  38. # the two calculations might be artificially close.
  39. raise PrecisionExhausted
  40. def _eigenvals_mpmath(M, multiple=False):
  41. """Compute eigenvalues using mpmath"""
  42. E, _ = _eigenvals_eigenvects_mpmath(M)
  43. result = [_sympify(x) for x in E]
  44. if multiple:
  45. return result
  46. return dict(Counter(result))
  47. def _eigenvects_mpmath(M):
  48. E, ER = _eigenvals_eigenvects_mpmath(M)
  49. result = []
  50. for i in range(M.rows):
  51. eigenval = _sympify(E[i])
  52. eigenvect = _sympify(ER[:, i])
  53. result.append((eigenval, 1, [eigenvect]))
  54. return result
  55. # This function is a candidate for caching if it gets implemented for matrices.
  56. def _eigenvals(
  57. M, error_when_incomplete=True, *, simplify=False, multiple=False,
  58. rational=False, **flags):
  59. r"""Compute eigenvalues of the matrix.
  60. Parameters
  61. ==========
  62. error_when_incomplete : bool, optional
  63. If it is set to ``True``, it will raise an error if not all
  64. eigenvalues are computed. This is caused by ``roots`` not returning
  65. a full list of eigenvalues.
  66. simplify : bool or function, optional
  67. If it is set to ``True``, it attempts to return the most
  68. simplified form of expressions returned by applying default
  69. simplification method in every routine.
  70. If it is set to ``False``, it will skip simplification in this
  71. particular routine to save computation resources.
  72. If a function is passed to, it will attempt to apply
  73. the particular function as simplification method.
  74. rational : bool, optional
  75. If it is set to ``True``, every floating point numbers would be
  76. replaced with rationals before computation. It can solve some
  77. issues of ``roots`` routine not working well with floats.
  78. multiple : bool, optional
  79. If it is set to ``True``, the result will be in the form of a
  80. list.
  81. If it is set to ``False``, the result will be in the form of a
  82. dictionary.
  83. Returns
  84. =======
  85. eigs : list or dict
  86. Eigenvalues of a matrix. The return format would be specified by
  87. the key ``multiple``.
  88. Raises
  89. ======
  90. MatrixError
  91. If not enough roots had got computed.
  92. NonSquareMatrixError
  93. If attempted to compute eigenvalues from a non-square matrix.
  94. Examples
  95. ========
  96. >>> from sympy import Matrix
  97. >>> M = Matrix(3, 3, [0, 1, 1, 1, 0, 0, 1, 1, 1])
  98. >>> M.eigenvals()
  99. {-1: 1, 0: 1, 2: 1}
  100. See Also
  101. ========
  102. MatrixDeterminant.charpoly
  103. eigenvects
  104. Notes
  105. =====
  106. Eigenvalues of a matrix $A$ can be computed by solving a matrix
  107. equation $\det(A - \lambda I) = 0$
  108. It's not always possible to return radical solutions for
  109. eigenvalues for matrices larger than $4, 4$ shape due to
  110. Abel-Ruffini theorem.
  111. If there is no radical solution is found for the eigenvalue,
  112. it may return eigenvalues in the form of
  113. :class:`sympy.polys.rootoftools.ComplexRootOf`.
  114. """
  115. if not M:
  116. if multiple:
  117. return []
  118. return {}
  119. if not M.is_square:
  120. raise NonSquareMatrixError("{} must be a square matrix.".format(M))
  121. if M._rep.domain not in (ZZ, QQ):
  122. # Skip this check for ZZ/QQ because it can be slow
  123. if all(x.is_number for x in M) and M.has(Float):
  124. return _eigenvals_mpmath(M, multiple=multiple)
  125. if rational:
  126. M = M.applyfunc(
  127. lambda x: nsimplify(x, rational=True) if x.has(Float) else x)
  128. if multiple:
  129. return _eigenvals_list(
  130. M, error_when_incomplete=error_when_incomplete, simplify=simplify,
  131. **flags)
  132. return _eigenvals_dict(
  133. M, error_when_incomplete=error_when_incomplete, simplify=simplify,
  134. **flags)
  135. eigenvals_error_message = \
  136. "It is not always possible to express the eigenvalues of a matrix " + \
  137. "of size 5x5 or higher in radicals. " + \
  138. "We have CRootOf, but domains other than the rationals are not " + \
  139. "currently supported. " + \
  140. "If there are no symbols in the matrix, " + \
  141. "it should still be possible to compute numeric approximations " + \
  142. "of the eigenvalues using " + \
  143. "M.evalf().eigenvals() or M.charpoly().nroots()."
  144. def _eigenvals_list(
  145. M, error_when_incomplete=True, simplify=False, **flags):
  146. iblocks = M.strongly_connected_components()
  147. all_eigs = []
  148. is_dom = M._rep.domain in (ZZ, QQ)
  149. for b in iblocks:
  150. # Fast path for a 1x1 block:
  151. if is_dom and len(b) == 1:
  152. index = b[0]
  153. val = M[index, index]
  154. all_eigs.append(val)
  155. continue
  156. block = M[b, b]
  157. if isinstance(simplify, FunctionType):
  158. charpoly = block.charpoly(simplify=simplify)
  159. else:
  160. charpoly = block.charpoly()
  161. eigs = roots(charpoly, multiple=True, **flags)
  162. if len(eigs) != block.rows:
  163. degree = int(charpoly.degree())
  164. f = charpoly.as_expr()
  165. x = charpoly.gen
  166. try:
  167. eigs = [CRootOf(f, x, idx) for idx in range(degree)]
  168. except NotImplementedError:
  169. if error_when_incomplete:
  170. raise MatrixError(eigenvals_error_message)
  171. else:
  172. eigs = []
  173. all_eigs += eigs
  174. if not simplify:
  175. return all_eigs
  176. if not isinstance(simplify, FunctionType):
  177. simplify = _simplify
  178. return [simplify(value) for value in all_eigs]
  179. def _eigenvals_dict(
  180. M, error_when_incomplete=True, simplify=False, **flags):
  181. iblocks = M.strongly_connected_components()
  182. all_eigs = {}
  183. is_dom = M._rep.domain in (ZZ, QQ)
  184. for b in iblocks:
  185. # Fast path for a 1x1 block:
  186. if is_dom and len(b) == 1:
  187. index = b[0]
  188. val = M[index, index]
  189. all_eigs[val] = all_eigs.get(val, 0) + 1
  190. continue
  191. block = M[b, b]
  192. if isinstance(simplify, FunctionType):
  193. charpoly = block.charpoly(simplify=simplify)
  194. else:
  195. charpoly = block.charpoly()
  196. eigs = roots(charpoly, multiple=False, **flags)
  197. if sum(eigs.values()) != block.rows:
  198. degree = int(charpoly.degree())
  199. f = charpoly.as_expr()
  200. x = charpoly.gen
  201. try:
  202. eigs = {CRootOf(f, x, idx): 1 for idx in range(degree)}
  203. except NotImplementedError:
  204. if error_when_incomplete:
  205. raise MatrixError(eigenvals_error_message)
  206. else:
  207. eigs = {}
  208. for k, v in eigs.items():
  209. if k in all_eigs:
  210. all_eigs[k] += v
  211. else:
  212. all_eigs[k] = v
  213. if not simplify:
  214. return all_eigs
  215. if not isinstance(simplify, FunctionType):
  216. simplify = _simplify
  217. return {simplify(key): value for key, value in all_eigs.items()}
  218. def _eigenspace(M, eigenval, iszerofunc=_iszero, simplify=False):
  219. """Get a basis for the eigenspace for a particular eigenvalue"""
  220. m = M - M.eye(M.rows) * eigenval
  221. ret = m.nullspace(iszerofunc=iszerofunc)
  222. # The nullspace for a real eigenvalue should be non-trivial.
  223. # If we didn't find an eigenvector, try once more a little harder
  224. if len(ret) == 0 and simplify:
  225. ret = m.nullspace(iszerofunc=iszerofunc, simplify=True)
  226. if len(ret) == 0:
  227. raise NotImplementedError(
  228. "Can't evaluate eigenvector for eigenvalue {}".format(eigenval))
  229. return ret
  230. def _eigenvects_DOM(M, **kwargs):
  231. DOM = DomainMatrix.from_Matrix(M, field=True, extension=True)
  232. DOM = DOM.to_dense()
  233. if DOM.domain != EX:
  234. rational, algebraic = dom_eigenvects(DOM)
  235. eigenvects = dom_eigenvects_to_sympy(
  236. rational, algebraic, M.__class__, **kwargs)
  237. eigenvects = sorted(eigenvects, key=lambda x: default_sort_key(x[0]))
  238. return eigenvects
  239. return None
  240. def _eigenvects_sympy(M, iszerofunc, simplify=True, **flags):
  241. eigenvals = M.eigenvals(rational=False, **flags)
  242. # Make sure that we have all roots in radical form
  243. for x in eigenvals:
  244. if x.has(CRootOf):
  245. raise MatrixError(
  246. "Eigenvector computation is not implemented if the matrix have "
  247. "eigenvalues in CRootOf form")
  248. eigenvals = sorted(eigenvals.items(), key=default_sort_key)
  249. ret = []
  250. for val, mult in eigenvals:
  251. vects = _eigenspace(M, val, iszerofunc=iszerofunc, simplify=simplify)
  252. ret.append((val, mult, vects))
  253. return ret
  254. # This functions is a candidate for caching if it gets implemented for matrices.
  255. def _eigenvects(M, error_when_incomplete=True, iszerofunc=_iszero, *, chop=False, **flags):
  256. """Compute eigenvectors of the matrix.
  257. Parameters
  258. ==========
  259. error_when_incomplete : bool, optional
  260. Raise an error when not all eigenvalues are computed. This is
  261. caused by ``roots`` not returning a full list of eigenvalues.
  262. iszerofunc : function, optional
  263. Specifies a zero testing function to be used in ``rref``.
  264. Default value is ``_iszero``, which uses SymPy's naive and fast
  265. default assumption handler.
  266. It can also accept any user-specified zero testing function, if it
  267. is formatted as a function which accepts a single symbolic argument
  268. and returns ``True`` if it is tested as zero and ``False`` if it
  269. is tested as non-zero, and ``None`` if it is undecidable.
  270. simplify : bool or function, optional
  271. If ``True``, ``as_content_primitive()`` will be used to tidy up
  272. normalization artifacts.
  273. It will also be used by the ``nullspace`` routine.
  274. chop : bool or positive number, optional
  275. If the matrix contains any Floats, they will be changed to Rationals
  276. for computation purposes, but the answers will be returned after
  277. being evaluated with evalf. The ``chop`` flag is passed to ``evalf``.
  278. When ``chop=True`` a default precision will be used; a number will
  279. be interpreted as the desired level of precision.
  280. Returns
  281. =======
  282. ret : [(eigenval, multiplicity, eigenspace), ...]
  283. A ragged list containing tuples of data obtained by ``eigenvals``
  284. and ``nullspace``.
  285. ``eigenspace`` is a list containing the ``eigenvector`` for each
  286. eigenvalue.
  287. ``eigenvector`` is a vector in the form of a ``Matrix``. e.g.
  288. a vector of length 3 is returned as ``Matrix([a_1, a_2, a_3])``.
  289. Raises
  290. ======
  291. NotImplementedError
  292. If failed to compute nullspace.
  293. Examples
  294. ========
  295. >>> from sympy import Matrix
  296. >>> M = Matrix(3, 3, [0, 1, 1, 1, 0, 0, 1, 1, 1])
  297. >>> M.eigenvects()
  298. [(-1, 1, [Matrix([
  299. [-1],
  300. [ 1],
  301. [ 0]])]), (0, 1, [Matrix([
  302. [ 0],
  303. [-1],
  304. [ 1]])]), (2, 1, [Matrix([
  305. [2/3],
  306. [1/3],
  307. [ 1]])])]
  308. See Also
  309. ========
  310. eigenvals
  311. MatrixSubspaces.nullspace
  312. """
  313. simplify = flags.get('simplify', True)
  314. primitive = flags.get('simplify', False)
  315. flags.pop('simplify', None) # remove this if it's there
  316. flags.pop('multiple', None) # remove this if it's there
  317. if not isinstance(simplify, FunctionType):
  318. simpfunc = _simplify if simplify else lambda x: x
  319. has_floats = M.has(Float)
  320. if has_floats:
  321. if all(x.is_number for x in M):
  322. return _eigenvects_mpmath(M)
  323. M = M.applyfunc(lambda x: nsimplify(x, rational=True))
  324. ret = _eigenvects_DOM(M)
  325. if ret is None:
  326. ret = _eigenvects_sympy(M, iszerofunc, simplify=simplify, **flags)
  327. if primitive:
  328. # if the primitive flag is set, get rid of any common
  329. # integer denominators
  330. def denom_clean(l):
  331. from sympy.polys.polytools import gcd
  332. return [(v / gcd(list(v))).applyfunc(simpfunc) for v in l]
  333. ret = [(val, mult, denom_clean(es)) for val, mult, es in ret]
  334. if has_floats:
  335. # if we had floats to start with, turn the eigenvectors to floats
  336. ret = [(val.evalf(chop=chop), mult, [v.evalf(chop=chop) for v in es])
  337. for val, mult, es in ret]
  338. return ret
  339. def _is_diagonalizable_with_eigen(M, reals_only=False):
  340. """See _is_diagonalizable. This function returns the bool along with the
  341. eigenvectors to avoid calculating them again in functions like
  342. ``diagonalize``."""
  343. if not M.is_square:
  344. return False, []
  345. eigenvecs = M.eigenvects(simplify=True)
  346. for val, mult, basis in eigenvecs:
  347. if reals_only and not val.is_real: # if we have a complex eigenvalue
  348. return False, eigenvecs
  349. if mult != len(basis): # if the geometric multiplicity doesn't equal the algebraic
  350. return False, eigenvecs
  351. return True, eigenvecs
  352. def _is_diagonalizable(M, reals_only=False, **kwargs):
  353. """Returns ``True`` if a matrix is diagonalizable.
  354. Parameters
  355. ==========
  356. reals_only : bool, optional
  357. If ``True``, it tests whether the matrix can be diagonalized
  358. to contain only real numbers on the diagonal.
  359. If ``False``, it tests whether the matrix can be diagonalized
  360. at all, even with numbers that may not be real.
  361. Examples
  362. ========
  363. Example of a diagonalizable matrix:
  364. >>> from sympy import Matrix
  365. >>> M = Matrix([[1, 2, 0], [0, 3, 0], [2, -4, 2]])
  366. >>> M.is_diagonalizable()
  367. True
  368. Example of a non-diagonalizable matrix:
  369. >>> M = Matrix([[0, 1], [0, 0]])
  370. >>> M.is_diagonalizable()
  371. False
  372. Example of a matrix that is diagonalized in terms of non-real entries:
  373. >>> M = Matrix([[0, 1], [-1, 0]])
  374. >>> M.is_diagonalizable(reals_only=False)
  375. True
  376. >>> M.is_diagonalizable(reals_only=True)
  377. False
  378. See Also
  379. ========
  380. is_diagonal
  381. diagonalize
  382. """
  383. if 'clear_cache' in kwargs:
  384. sympy_deprecation_warning(
  385. """
  386. The 'clear_cache' keyword to Matrix.is_diagonalizable is
  387. deprecated. It does nothing and should be omitted.
  388. """,
  389. deprecated_since_version="1.4",
  390. active_deprecations_target="deprecated-matrix-is_diagonalizable-cache",
  391. stacklevel=4,
  392. )
  393. if 'clear_subproducts' in kwargs:
  394. sympy_deprecation_warning(
  395. """
  396. The 'clear_subproducts' keyword to Matrix.is_diagonalizable is
  397. deprecated. It does nothing and should be omitted.
  398. """,
  399. deprecated_since_version="1.4",
  400. active_deprecations_target="deprecated-matrix-is_diagonalizable-cache",
  401. stacklevel=4,
  402. )
  403. if not M.is_square:
  404. return False
  405. if all(e.is_real for e in M) and M.is_symmetric():
  406. return True
  407. if all(e.is_complex for e in M) and M.is_hermitian:
  408. return True
  409. return _is_diagonalizable_with_eigen(M, reals_only=reals_only)[0]
  410. #G&VL, Matrix Computations, Algo 5.4.2
  411. def _householder_vector(x):
  412. if not x.cols == 1:
  413. raise ValueError("Input must be a column matrix")
  414. v = x.copy()
  415. v_plus = x.copy()
  416. v_minus = x.copy()
  417. q = x[0, 0] / abs(x[0, 0])
  418. norm_x = x.norm()
  419. v_plus[0, 0] = x[0, 0] + q * norm_x
  420. v_minus[0, 0] = x[0, 0] - q * norm_x
  421. if x[1:, 0].norm() == 0:
  422. bet = 0
  423. v[0, 0] = 1
  424. else:
  425. if v_plus.norm() <= v_minus.norm():
  426. v = v_plus
  427. else:
  428. v = v_minus
  429. v = v / v[0]
  430. bet = 2 / (v.norm() ** 2)
  431. return v, bet
  432. def _bidiagonal_decmp_hholder(M):
  433. m = M.rows
  434. n = M.cols
  435. A = M.as_mutable()
  436. U, V = A.eye(m), A.eye(n)
  437. for i in range(min(m, n)):
  438. v, bet = _householder_vector(A[i:, i])
  439. hh_mat = A.eye(m - i) - bet * v * v.H
  440. A[i:, i:] = hh_mat * A[i:, i:]
  441. temp = A.eye(m)
  442. temp[i:, i:] = hh_mat
  443. U = U * temp
  444. if i + 1 <= n - 2:
  445. v, bet = _householder_vector(A[i, i+1:].T)
  446. hh_mat = A.eye(n - i - 1) - bet * v * v.H
  447. A[i:, i+1:] = A[i:, i+1:] * hh_mat
  448. temp = A.eye(n)
  449. temp[i+1:, i+1:] = hh_mat
  450. V = temp * V
  451. return U, A, V
  452. def _eval_bidiag_hholder(M):
  453. m = M.rows
  454. n = M.cols
  455. A = M.as_mutable()
  456. for i in range(min(m, n)):
  457. v, bet = _householder_vector(A[i:, i])
  458. hh_mat = A.eye(m-i) - bet * v * v.H
  459. A[i:, i:] = hh_mat * A[i:, i:]
  460. if i + 1 <= n - 2:
  461. v, bet = _householder_vector(A[i, i+1:].T)
  462. hh_mat = A.eye(n - i - 1) - bet * v * v.H
  463. A[i:, i+1:] = A[i:, i+1:] * hh_mat
  464. return A
  465. def _bidiagonal_decomposition(M, upper=True):
  466. """
  467. Returns $(U,B,V.H)$ for
  468. $$A = UBV^{H}$$
  469. where $A$ is the input matrix, and $B$ is its Bidiagonalized form
  470. Note: Bidiagonal Computation can hang for symbolic matrices.
  471. Parameters
  472. ==========
  473. upper : bool. Whether to do upper bidiagnalization or lower.
  474. True for upper and False for lower.
  475. References
  476. ==========
  477. .. [1] Algorithm 5.4.2, Matrix computations by Golub and Van Loan, 4th edition
  478. .. [2] Complex Matrix Bidiagonalization, https://github.com/vslobody/Householder-Bidiagonalization
  479. """
  480. if not isinstance(upper, bool):
  481. raise ValueError("upper must be a boolean")
  482. if upper:
  483. return _bidiagonal_decmp_hholder(M)
  484. X = _bidiagonal_decmp_hholder(M.H)
  485. return X[2].H, X[1].H, X[0].H
  486. def _bidiagonalize(M, upper=True):
  487. """
  488. Returns $B$, the Bidiagonalized form of the input matrix.
  489. Note: Bidiagonal Computation can hang for symbolic matrices.
  490. Parameters
  491. ==========
  492. upper : bool. Whether to do upper bidiagnalization or lower.
  493. True for upper and False for lower.
  494. References
  495. ==========
  496. .. [1] Algorithm 5.4.2, Matrix computations by Golub and Van Loan, 4th edition
  497. .. [2] Complex Matrix Bidiagonalization : https://github.com/vslobody/Householder-Bidiagonalization
  498. """
  499. if not isinstance(upper, bool):
  500. raise ValueError("upper must be a boolean")
  501. if upper:
  502. return _eval_bidiag_hholder(M)
  503. return _eval_bidiag_hholder(M.H).H
  504. def _diagonalize(M, reals_only=False, sort=False, normalize=False):
  505. """
  506. Return (P, D), where D is diagonal and
  507. D = P^-1 * M * P
  508. where M is current matrix.
  509. Parameters
  510. ==========
  511. reals_only : bool. Whether to throw an error if complex numbers are need
  512. to diagonalize. (Default: False)
  513. sort : bool. Sort the eigenvalues along the diagonal. (Default: False)
  514. normalize : bool. If True, normalize the columns of P. (Default: False)
  515. Examples
  516. ========
  517. >>> from sympy import Matrix
  518. >>> M = Matrix(3, 3, [1, 2, 0, 0, 3, 0, 2, -4, 2])
  519. >>> M
  520. Matrix([
  521. [1, 2, 0],
  522. [0, 3, 0],
  523. [2, -4, 2]])
  524. >>> (P, D) = M.diagonalize()
  525. >>> D
  526. Matrix([
  527. [1, 0, 0],
  528. [0, 2, 0],
  529. [0, 0, 3]])
  530. >>> P
  531. Matrix([
  532. [-1, 0, -1],
  533. [ 0, 0, -1],
  534. [ 2, 1, 2]])
  535. >>> P.inv() * M * P
  536. Matrix([
  537. [1, 0, 0],
  538. [0, 2, 0],
  539. [0, 0, 3]])
  540. See Also
  541. ========
  542. is_diagonal
  543. is_diagonalizable
  544. """
  545. if not M.is_square:
  546. raise NonSquareMatrixError()
  547. is_diagonalizable, eigenvecs = _is_diagonalizable_with_eigen(M,
  548. reals_only=reals_only)
  549. if not is_diagonalizable:
  550. raise MatrixError("Matrix is not diagonalizable")
  551. if sort:
  552. eigenvecs = sorted(eigenvecs, key=default_sort_key)
  553. p_cols, diag = [], []
  554. for val, mult, basis in eigenvecs:
  555. diag += [val] * mult
  556. p_cols += basis
  557. if normalize:
  558. p_cols = [v / v.norm() for v in p_cols]
  559. return M.hstack(*p_cols), M.diag(*diag)
  560. def _fuzzy_positive_definite(M):
  561. positive_diagonals = M._has_positive_diagonals()
  562. if positive_diagonals is False:
  563. return False
  564. if positive_diagonals and M.is_strongly_diagonally_dominant:
  565. return True
  566. return None
  567. def _fuzzy_positive_semidefinite(M):
  568. nonnegative_diagonals = M._has_nonnegative_diagonals()
  569. if nonnegative_diagonals is False:
  570. return False
  571. if nonnegative_diagonals and M.is_weakly_diagonally_dominant:
  572. return True
  573. return None
  574. def _is_positive_definite(M):
  575. if not M.is_hermitian:
  576. if not M.is_square:
  577. return False
  578. M = M + M.H
  579. fuzzy = _fuzzy_positive_definite(M)
  580. if fuzzy is not None:
  581. return fuzzy
  582. return _is_positive_definite_GE(M)
  583. def _is_positive_semidefinite(M):
  584. if not M.is_hermitian:
  585. if not M.is_square:
  586. return False
  587. M = M + M.H
  588. fuzzy = _fuzzy_positive_semidefinite(M)
  589. if fuzzy is not None:
  590. return fuzzy
  591. return _is_positive_semidefinite_cholesky(M)
  592. def _is_negative_definite(M):
  593. return _is_positive_definite(-M)
  594. def _is_negative_semidefinite(M):
  595. return _is_positive_semidefinite(-M)
  596. def _is_indefinite(M):
  597. if M.is_hermitian:
  598. eigen = M.eigenvals()
  599. args1 = [x.is_positive for x in eigen.keys()]
  600. any_positive = fuzzy_or(args1)
  601. args2 = [x.is_negative for x in eigen.keys()]
  602. any_negative = fuzzy_or(args2)
  603. return fuzzy_and([any_positive, any_negative])
  604. elif M.is_square:
  605. return (M + M.H).is_indefinite
  606. return False
  607. def _is_positive_definite_GE(M):
  608. """A division-free gaussian elimination method for testing
  609. positive-definiteness."""
  610. M = M.as_mutable()
  611. size = M.rows
  612. for i in range(size):
  613. is_positive = M[i, i].is_positive
  614. if is_positive is not True:
  615. return is_positive
  616. for j in range(i+1, size):
  617. M[j, i+1:] = M[i, i] * M[j, i+1:] - M[j, i] * M[i, i+1:]
  618. return True
  619. def _is_positive_semidefinite_cholesky(M):
  620. """Uses Cholesky factorization with complete pivoting
  621. References
  622. ==========
  623. .. [1] http://eprints.ma.man.ac.uk/1199/1/covered/MIMS_ep2008_116.pdf
  624. .. [2] https://www.value-at-risk.net/cholesky-factorization/
  625. """
  626. M = M.as_mutable()
  627. for k in range(M.rows):
  628. diags = [M[i, i] for i in range(k, M.rows)]
  629. pivot, pivot_val, nonzero, _ = _find_reasonable_pivot(diags)
  630. if nonzero:
  631. return None
  632. if pivot is None:
  633. for i in range(k+1, M.rows):
  634. for j in range(k, M.cols):
  635. iszero = M[i, j].is_zero
  636. if iszero is None:
  637. return None
  638. elif iszero is False:
  639. return False
  640. return True
  641. if M[k, k].is_negative or pivot_val.is_negative:
  642. return False
  643. elif not (M[k, k].is_nonnegative and pivot_val.is_nonnegative):
  644. return None
  645. if pivot > 0:
  646. M.col_swap(k, k+pivot)
  647. M.row_swap(k, k+pivot)
  648. M[k, k] = sqrt(M[k, k])
  649. M[k, k+1:] /= M[k, k]
  650. M[k+1:, k+1:] -= M[k, k+1:].H * M[k, k+1:]
  651. return M[-1, -1].is_nonnegative
  652. _doc_positive_definite = \
  653. r"""Finds out the definiteness of a matrix.
  654. Explanation
  655. ===========
  656. A square real matrix $A$ is:
  657. - A positive definite matrix if $x^T A x > 0$
  658. for all non-zero real vectors $x$.
  659. - A positive semidefinite matrix if $x^T A x \geq 0$
  660. for all non-zero real vectors $x$.
  661. - A negative definite matrix if $x^T A x < 0$
  662. for all non-zero real vectors $x$.
  663. - A negative semidefinite matrix if $x^T A x \leq 0$
  664. for all non-zero real vectors $x$.
  665. - An indefinite matrix if there exists non-zero real vectors
  666. $x, y$ with $x^T A x > 0 > y^T A y$.
  667. A square complex matrix $A$ is:
  668. - A positive definite matrix if $\text{re}(x^H A x) > 0$
  669. for all non-zero complex vectors $x$.
  670. - A positive semidefinite matrix if $\text{re}(x^H A x) \geq 0$
  671. for all non-zero complex vectors $x$.
  672. - A negative definite matrix if $\text{re}(x^H A x) < 0$
  673. for all non-zero complex vectors $x$.
  674. - A negative semidefinite matrix if $\text{re}(x^H A x) \leq 0$
  675. for all non-zero complex vectors $x$.
  676. - An indefinite matrix if there exists non-zero complex vectors
  677. $x, y$ with $\text{re}(x^H A x) > 0 > \text{re}(y^H A y)$.
  678. A matrix need not be symmetric or hermitian to be positive definite.
  679. - A real non-symmetric matrix is positive definite if and only if
  680. $\frac{A + A^T}{2}$ is positive definite.
  681. - A complex non-hermitian matrix is positive definite if and only if
  682. $\frac{A + A^H}{2}$ is positive definite.
  683. And this extension can apply for all the definitions above.
  684. However, for complex cases, you can restrict the definition of
  685. $\text{re}(x^H A x) > 0$ to $x^H A x > 0$ and require the matrix
  686. to be hermitian.
  687. But we do not present this restriction for computation because you
  688. can check ``M.is_hermitian`` independently with this and use
  689. the same procedure.
  690. Examples
  691. ========
  692. An example of symmetric positive definite matrix:
  693. .. plot::
  694. :context: reset
  695. :format: doctest
  696. :include-source: True
  697. >>> from sympy import Matrix, symbols
  698. >>> from sympy.plotting import plot3d
  699. >>> a, b = symbols('a b')
  700. >>> x = Matrix([a, b])
  701. >>> A = Matrix([[1, 0], [0, 1]])
  702. >>> A.is_positive_definite
  703. True
  704. >>> A.is_positive_semidefinite
  705. True
  706. >>> p = plot3d((x.T*A*x)[0, 0], (a, -1, 1), (b, -1, 1))
  707. An example of symmetric positive semidefinite matrix:
  708. .. plot::
  709. :context: close-figs
  710. :format: doctest
  711. :include-source: True
  712. >>> A = Matrix([[1, -1], [-1, 1]])
  713. >>> A.is_positive_definite
  714. False
  715. >>> A.is_positive_semidefinite
  716. True
  717. >>> p = plot3d((x.T*A*x)[0, 0], (a, -1, 1), (b, -1, 1))
  718. An example of symmetric negative definite matrix:
  719. .. plot::
  720. :context: close-figs
  721. :format: doctest
  722. :include-source: True
  723. >>> A = Matrix([[-1, 0], [0, -1]])
  724. >>> A.is_negative_definite
  725. True
  726. >>> A.is_negative_semidefinite
  727. True
  728. >>> A.is_indefinite
  729. False
  730. >>> p = plot3d((x.T*A*x)[0, 0], (a, -1, 1), (b, -1, 1))
  731. An example of symmetric indefinite matrix:
  732. .. plot::
  733. :context: close-figs
  734. :format: doctest
  735. :include-source: True
  736. >>> A = Matrix([[1, 2], [2, -1]])
  737. >>> A.is_indefinite
  738. True
  739. >>> p = plot3d((x.T*A*x)[0, 0], (a, -1, 1), (b, -1, 1))
  740. An example of non-symmetric positive definite matrix.
  741. .. plot::
  742. :context: close-figs
  743. :format: doctest
  744. :include-source: True
  745. >>> A = Matrix([[1, 2], [-2, 1]])
  746. >>> A.is_positive_definite
  747. True
  748. >>> A.is_positive_semidefinite
  749. True
  750. >>> p = plot3d((x.T*A*x)[0, 0], (a, -1, 1), (b, -1, 1))
  751. Notes
  752. =====
  753. Although some people trivialize the definition of positive definite
  754. matrices only for symmetric or hermitian matrices, this restriction
  755. is not correct because it does not classify all instances of
  756. positive definite matrices from the definition $x^T A x > 0$ or
  757. $\text{re}(x^H A x) > 0$.
  758. For instance, ``Matrix([[1, 2], [-2, 1]])`` presented in
  759. the example above is an example of real positive definite matrix
  760. that is not symmetric.
  761. However, since the following formula holds true;
  762. .. math::
  763. \text{re}(x^H A x) > 0 \iff
  764. \text{re}(x^H \frac{A + A^H}{2} x) > 0
  765. We can classify all positive definite matrices that may or may not
  766. be symmetric or hermitian by transforming the matrix to
  767. $\frac{A + A^T}{2}$ or $\frac{A + A^H}{2}$
  768. (which is guaranteed to be always real symmetric or complex
  769. hermitian) and we can defer most of the studies to symmetric or
  770. hermitian positive definite matrices.
  771. But it is a different problem for the existance of Cholesky
  772. decomposition. Because even though a non symmetric or a non
  773. hermitian matrix can be positive definite, Cholesky or LDL
  774. decomposition does not exist because the decompositions require the
  775. matrix to be symmetric or hermitian.
  776. References
  777. ==========
  778. .. [1] https://en.wikipedia.org/wiki/Definiteness_of_a_matrix#Eigenvalues
  779. .. [2] http://mathworld.wolfram.com/PositiveDefiniteMatrix.html
  780. .. [3] Johnson, C. R. "Positive Definite Matrices." Amer.
  781. Math. Monthly 77, 259-264 1970.
  782. """
  783. _is_positive_definite.__doc__ = _doc_positive_definite
  784. _is_positive_semidefinite.__doc__ = _doc_positive_definite
  785. _is_negative_definite.__doc__ = _doc_positive_definite
  786. _is_negative_semidefinite.__doc__ = _doc_positive_definite
  787. _is_indefinite.__doc__ = _doc_positive_definite
  788. def _jordan_form(M, calc_transform=True, *, chop=False):
  789. """Return $(P, J)$ where $J$ is a Jordan block
  790. matrix and $P$ is a matrix such that $M = P J P^{-1}$
  791. Parameters
  792. ==========
  793. calc_transform : bool
  794. If ``False``, then only $J$ is returned.
  795. chop : bool
  796. All matrices are converted to exact types when computing
  797. eigenvalues and eigenvectors. As a result, there may be
  798. approximation errors. If ``chop==True``, these errors
  799. will be truncated.
  800. Examples
  801. ========
  802. >>> from sympy import Matrix
  803. >>> M = Matrix([[ 6, 5, -2, -3], [-3, -1, 3, 3], [ 2, 1, -2, -3], [-1, 1, 5, 5]])
  804. >>> P, J = M.jordan_form()
  805. >>> J
  806. Matrix([
  807. [2, 1, 0, 0],
  808. [0, 2, 0, 0],
  809. [0, 0, 2, 1],
  810. [0, 0, 0, 2]])
  811. See Also
  812. ========
  813. jordan_block
  814. """
  815. if not M.is_square:
  816. raise NonSquareMatrixError("Only square matrices have Jordan forms")
  817. mat = M
  818. has_floats = M.has(Float)
  819. if has_floats:
  820. try:
  821. max_prec = max(term._prec for term in M.values() if isinstance(term, Float))
  822. except ValueError:
  823. # if no term in the matrix is explicitly a Float calling max()
  824. # will throw a error so setting max_prec to default value of 53
  825. max_prec = 53
  826. # setting minimum max_dps to 15 to prevent loss of precision in
  827. # matrix containing non evaluated expressions
  828. max_dps = max(prec_to_dps(max_prec), 15)
  829. def restore_floats(*args):
  830. """If ``has_floats`` is `True`, cast all ``args`` as
  831. matrices of floats."""
  832. if has_floats:
  833. args = [m.evalf(n=max_dps, chop=chop) for m in args]
  834. if len(args) == 1:
  835. return args[0]
  836. return args
  837. # cache calculations for some speedup
  838. mat_cache = {}
  839. def eig_mat(val, pow):
  840. """Cache computations of ``(M - val*I)**pow`` for quick
  841. retrieval"""
  842. if (val, pow) in mat_cache:
  843. return mat_cache[(val, pow)]
  844. if (val, pow - 1) in mat_cache:
  845. mat_cache[(val, pow)] = mat_cache[(val, pow - 1)].multiply(
  846. mat_cache[(val, 1)], dotprodsimp=None)
  847. else:
  848. mat_cache[(val, pow)] = (mat - val*M.eye(M.rows)).pow(pow)
  849. return mat_cache[(val, pow)]
  850. # helper functions
  851. def nullity_chain(val, algebraic_multiplicity):
  852. """Calculate the sequence [0, nullity(E), nullity(E**2), ...]
  853. until it is constant where ``E = M - val*I``"""
  854. # mat.rank() is faster than computing the null space,
  855. # so use the rank-nullity theorem
  856. cols = M.cols
  857. ret = [0]
  858. nullity = cols - eig_mat(val, 1).rank()
  859. i = 2
  860. while nullity != ret[-1]:
  861. ret.append(nullity)
  862. if nullity == algebraic_multiplicity:
  863. break
  864. nullity = cols - eig_mat(val, i).rank()
  865. i += 1
  866. # Due to issues like #7146 and #15872, SymPy sometimes
  867. # gives the wrong rank. In this case, raise an error
  868. # instead of returning an incorrect matrix
  869. if nullity < ret[-1] or nullity > algebraic_multiplicity:
  870. raise MatrixError(
  871. "SymPy had encountered an inconsistent "
  872. "result while computing Jordan block: "
  873. "{}".format(M))
  874. return ret
  875. def blocks_from_nullity_chain(d):
  876. """Return a list of the size of each Jordan block.
  877. If d_n is the nullity of E**n, then the number
  878. of Jordan blocks of size n is
  879. 2*d_n - d_(n-1) - d_(n+1)"""
  880. # d[0] is always the number of columns, so skip past it
  881. mid = [2*d[n] - d[n - 1] - d[n + 1] for n in range(1, len(d) - 1)]
  882. # d is assumed to plateau with "d[ len(d) ] == d[-1]", so
  883. # 2*d_n - d_(n-1) - d_(n+1) == d_n - d_(n-1)
  884. end = [d[-1] - d[-2]] if len(d) > 1 else [d[0]]
  885. return mid + end
  886. def pick_vec(small_basis, big_basis):
  887. """Picks a vector from big_basis that isn't in
  888. the subspace spanned by small_basis"""
  889. if len(small_basis) == 0:
  890. return big_basis[0]
  891. for v in big_basis:
  892. _, pivots = M.hstack(*(small_basis + [v])).echelon_form(
  893. with_pivots=True)
  894. if pivots[-1] == len(small_basis):
  895. return v
  896. # roots doesn't like Floats, so replace them with Rationals
  897. if has_floats:
  898. mat = mat.applyfunc(lambda x: nsimplify(x, rational=True))
  899. # first calculate the jordan block structure
  900. eigs = mat.eigenvals()
  901. # Make sure that we have all roots in radical form
  902. for x in eigs:
  903. if x.has(CRootOf):
  904. raise MatrixError(
  905. "Jordan normal form is not implemented if the matrix have "
  906. "eigenvalues in CRootOf form")
  907. # most matrices have distinct eigenvalues
  908. # and so are diagonalizable. In this case, don't
  909. # do extra work!
  910. if len(eigs.keys()) == mat.cols:
  911. blocks = list(sorted(eigs.keys(), key=default_sort_key))
  912. jordan_mat = mat.diag(*blocks)
  913. if not calc_transform:
  914. return restore_floats(jordan_mat)
  915. jordan_basis = [eig_mat(eig, 1).nullspace()[0]
  916. for eig in blocks]
  917. basis_mat = mat.hstack(*jordan_basis)
  918. return restore_floats(basis_mat, jordan_mat)
  919. block_structure = []
  920. for eig in sorted(eigs.keys(), key=default_sort_key):
  921. algebraic_multiplicity = eigs[eig]
  922. chain = nullity_chain(eig, algebraic_multiplicity)
  923. block_sizes = blocks_from_nullity_chain(chain)
  924. # if block_sizes = = [a, b, c, ...], then the number of
  925. # Jordan blocks of size 1 is a, of size 2 is b, etc.
  926. # create an array that has (eig, block_size) with one
  927. # entry for each block
  928. size_nums = [(i+1, num) for i, num in enumerate(block_sizes)]
  929. # we expect larger Jordan blocks to come earlier
  930. size_nums.reverse()
  931. block_structure.extend(
  932. [(eig, size) for size, num in size_nums for _ in range(num)])
  933. jordan_form_size = sum(size for eig, size in block_structure)
  934. if jordan_form_size != M.rows:
  935. raise MatrixError(
  936. "SymPy had encountered an inconsistent result while "
  937. "computing Jordan block. : {}".format(M))
  938. blocks = (mat.jordan_block(size=size, eigenvalue=eig) for eig, size in block_structure)
  939. jordan_mat = mat.diag(*blocks)
  940. if not calc_transform:
  941. return restore_floats(jordan_mat)
  942. # For each generalized eigenspace, calculate a basis.
  943. # We start by looking for a vector in null( (A - eig*I)**n )
  944. # which isn't in null( (A - eig*I)**(n-1) ) where n is
  945. # the size of the Jordan block
  946. #
  947. # Ideally we'd just loop through block_structure and
  948. # compute each generalized eigenspace. However, this
  949. # causes a lot of unneeded computation. Instead, we
  950. # go through the eigenvalues separately, since we know
  951. # their generalized eigenspaces must have bases that
  952. # are linearly independent.
  953. jordan_basis = []
  954. for eig in sorted(eigs.keys(), key=default_sort_key):
  955. eig_basis = []
  956. for block_eig, size in block_structure:
  957. if block_eig != eig:
  958. continue
  959. null_big = (eig_mat(eig, size)).nullspace()
  960. null_small = (eig_mat(eig, size - 1)).nullspace()
  961. # we want to pick something that is in the big basis
  962. # and not the small, but also something that is independent
  963. # of any other generalized eigenvectors from a different
  964. # generalized eigenspace sharing the same eigenvalue.
  965. vec = pick_vec(null_small + eig_basis, null_big)
  966. new_vecs = [eig_mat(eig, i).multiply(vec, dotprodsimp=None)
  967. for i in range(size)]
  968. eig_basis.extend(new_vecs)
  969. jordan_basis.extend(reversed(new_vecs))
  970. basis_mat = mat.hstack(*jordan_basis)
  971. return restore_floats(basis_mat, jordan_mat)
  972. def _left_eigenvects(M, **flags):
  973. """Returns left eigenvectors and eigenvalues.
  974. This function returns the list of triples (eigenval, multiplicity,
  975. basis) for the left eigenvectors. Options are the same as for
  976. eigenvects(), i.e. the ``**flags`` arguments gets passed directly to
  977. eigenvects().
  978. Examples
  979. ========
  980. >>> from sympy import Matrix
  981. >>> M = Matrix([[0, 1, 1], [1, 0, 0], [1, 1, 1]])
  982. >>> M.eigenvects()
  983. [(-1, 1, [Matrix([
  984. [-1],
  985. [ 1],
  986. [ 0]])]), (0, 1, [Matrix([
  987. [ 0],
  988. [-1],
  989. [ 1]])]), (2, 1, [Matrix([
  990. [2/3],
  991. [1/3],
  992. [ 1]])])]
  993. >>> M.left_eigenvects()
  994. [(-1, 1, [Matrix([[-2, 1, 1]])]), (0, 1, [Matrix([[-1, -1, 1]])]), (2,
  995. 1, [Matrix([[1, 1, 1]])])]
  996. """
  997. eigs = M.transpose().eigenvects(**flags)
  998. return [(val, mult, [l.transpose() for l in basis]) for val, mult, basis in eigs]
  999. def _singular_values(M):
  1000. """Compute the singular values of a Matrix
  1001. Examples
  1002. ========
  1003. >>> from sympy import Matrix, Symbol
  1004. >>> x = Symbol('x', real=True)
  1005. >>> M = Matrix([[0, 1, 0], [0, x, 0], [-1, 0, 0]])
  1006. >>> M.singular_values()
  1007. [sqrt(x**2 + 1), 1, 0]
  1008. See Also
  1009. ========
  1010. condition_number
  1011. """
  1012. if M.rows >= M.cols:
  1013. valmultpairs = M.H.multiply(M).eigenvals()
  1014. else:
  1015. valmultpairs = M.multiply(M.H).eigenvals()
  1016. # Expands result from eigenvals into a simple list
  1017. vals = []
  1018. for k, v in valmultpairs.items():
  1019. vals += [sqrt(k)] * v # dangerous! same k in several spots!
  1020. # Pad with zeros if singular values are computed in reverse way,
  1021. # to give consistent format.
  1022. if len(vals) < M.cols:
  1023. vals += [M.zero] * (M.cols - len(vals))
  1024. # sort them in descending order
  1025. vals.sort(reverse=True, key=default_sort_key)
  1026. return vals