_decomp.py 61 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651
  1. #
  2. # Author: Pearu Peterson, March 2002
  3. #
  4. # additions by Travis Oliphant, March 2002
  5. # additions by Eric Jones, June 2002
  6. # additions by Johannes Loehnert, June 2006
  7. # additions by Bart Vandereycken, June 2006
  8. # additions by Andrew D Straw, May 2007
  9. # additions by Tiziano Zito, November 2008
  10. #
  11. # April 2010: Functions for LU, QR, SVD, Schur, and Cholesky decompositions
  12. # were moved to their own files. Still in this file are functions for
  13. # eigenstuff and for the Hessenberg form.
  14. __all__ = ['eig', 'eigvals', 'eigh', 'eigvalsh',
  15. 'eig_banded', 'eigvals_banded',
  16. 'eigh_tridiagonal', 'eigvalsh_tridiagonal', 'hessenberg', 'cdf2rdf']
  17. import numpy as np
  18. from numpy import (array, isfinite, inexact, nonzero, iscomplexobj,
  19. flatnonzero, conj, asarray, argsort, empty,
  20. iscomplex, zeros, einsum, eye, inf)
  21. # Local imports
  22. from scipy._lib._util import _asarray_validated, _apply_over_batch
  23. from ._misc import LinAlgError, _datacopied, norm
  24. from .lapack import get_lapack_funcs, _compute_lwork
  25. _I = np.array(1j, dtype='F')
  26. def _make_complex_eigvecs(w, vin, dtype):
  27. """
  28. Produce complex-valued eigenvectors from LAPACK DGGEV real-valued output
  29. """
  30. # - see LAPACK man page DGGEV at ALPHAI
  31. v = np.array(vin, dtype=dtype)
  32. m = (w.imag > 0)
  33. m[:-1] |= (w.imag[1:] < 0) # workaround for LAPACK bug, cf. ticket #709
  34. for i in flatnonzero(m):
  35. v.imag[:, i] = vin[:, i+1]
  36. conj(v[:, i], v[:, i+1])
  37. return v
  38. def _make_eigvals(alpha, beta, homogeneous_eigvals):
  39. if homogeneous_eigvals:
  40. if beta is None:
  41. return np.vstack((alpha, np.ones_like(alpha)))
  42. else:
  43. return np.vstack((alpha, beta))
  44. else:
  45. if beta is None:
  46. return alpha
  47. else:
  48. w = np.empty_like(alpha)
  49. alpha_zero = (alpha == 0)
  50. beta_zero = (beta == 0)
  51. beta_nonzero = ~beta_zero
  52. w[beta_nonzero] = alpha[beta_nonzero]/beta[beta_nonzero]
  53. # Use np.inf for complex values too since
  54. # 1/np.inf = 0, i.e., it correctly behaves as projective
  55. # infinity.
  56. w[~alpha_zero & beta_zero] = np.inf
  57. if np.all(alpha.imag == 0):
  58. w[alpha_zero & beta_zero] = np.nan
  59. else:
  60. w[alpha_zero & beta_zero] = complex(np.nan, np.nan)
  61. return w
  62. def _geneig(a1, b1, left, right, overwrite_a, overwrite_b,
  63. homogeneous_eigvals):
  64. ggev, = get_lapack_funcs(('ggev',), (a1, b1))
  65. cvl, cvr = left, right
  66. res = ggev(a1, b1, lwork=-1)
  67. lwork = res[-2][0].real.astype(np.int_)
  68. if ggev.typecode in 'cz':
  69. alpha, beta, vl, vr, work, info = ggev(a1, b1, cvl, cvr, lwork,
  70. overwrite_a, overwrite_b)
  71. w = _make_eigvals(alpha, beta, homogeneous_eigvals)
  72. else:
  73. alphar, alphai, beta, vl, vr, work, info = ggev(a1, b1, cvl, cvr,
  74. lwork, overwrite_a,
  75. overwrite_b)
  76. alpha = alphar + _I * alphai
  77. w = _make_eigvals(alpha, beta, homogeneous_eigvals)
  78. _check_info(info, 'generalized eig algorithm (ggev)')
  79. only_real = np.all(w.imag == 0.0)
  80. if not (ggev.typecode in 'cz' or only_real):
  81. t = w.dtype.char
  82. if left:
  83. vl = _make_complex_eigvecs(w, vl, t)
  84. if right:
  85. vr = _make_complex_eigvecs(w, vr, t)
  86. # the eigenvectors returned by the lapack function are NOT normalized
  87. for i in range(vr.shape[0]):
  88. if right:
  89. vr[:, i] /= norm(vr[:, i])
  90. if left:
  91. vl[:, i] /= norm(vl[:, i])
  92. if not (left or right):
  93. return w
  94. if left:
  95. if right:
  96. return w, vl, vr
  97. return w, vl
  98. return w, vr
  99. @_apply_over_batch(('a', 2), ('b', 2))
  100. def eig(a, b=None, left=False, right=True, overwrite_a=False,
  101. overwrite_b=False, check_finite=True, homogeneous_eigvals=False):
  102. """
  103. Solve an ordinary or generalized eigenvalue problem of a square matrix.
  104. Find eigenvalues w and right or left eigenvectors of a general matrix::
  105. a vr[:,i] = w[i] b vr[:,i]
  106. a.H vl[:,i] = w[i].conj() b.H vl[:,i]
  107. where ``.H`` is the Hermitian conjugation.
  108. Parameters
  109. ----------
  110. a : (M, M) array_like
  111. A complex or real matrix whose eigenvalues and eigenvectors
  112. will be computed.
  113. b : (M, M) array_like, optional
  114. Right-hand side matrix in a generalized eigenvalue problem.
  115. Default is None, identity matrix is assumed.
  116. left : bool, optional
  117. Whether to calculate and return left eigenvectors. Default is False.
  118. right : bool, optional
  119. Whether to calculate and return right eigenvectors. Default is True.
  120. overwrite_a : bool, optional
  121. Whether to overwrite `a`; may improve performance. Default is False.
  122. overwrite_b : bool, optional
  123. Whether to overwrite `b`; may improve performance. Default is False.
  124. check_finite : bool, optional
  125. Whether to check that the input matrices contain only finite numbers.
  126. Disabling may give a performance gain, but may result in problems
  127. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  128. homogeneous_eigvals : bool, optional
  129. If True, return the eigenvalues in homogeneous coordinates.
  130. In this case ``w`` is a (2, M) array so that::
  131. w[1,i] a vr[:,i] = w[0,i] b vr[:,i]
  132. Default is False.
  133. Returns
  134. -------
  135. w : (M,) or (2, M) double or complex ndarray
  136. The eigenvalues, each repeated according to its
  137. multiplicity. The shape is (M,) unless
  138. ``homogeneous_eigvals=True``.
  139. vl : (M, M) double or complex ndarray
  140. The left eigenvector corresponding to the eigenvalue
  141. ``w[i]`` is the column ``vl[:,i]``. Only returned if ``left=True``.
  142. The left eigenvector is not normalized.
  143. vr : (M, M) double or complex ndarray
  144. The normalized right eigenvector corresponding to the eigenvalue
  145. ``w[i]`` is the column ``vr[:,i]``. Only returned if ``right=True``.
  146. Raises
  147. ------
  148. LinAlgError
  149. If eigenvalue computation does not converge.
  150. See Also
  151. --------
  152. eigvals : eigenvalues of general arrays
  153. eigh : Eigenvalues and right eigenvectors for symmetric/Hermitian arrays.
  154. eig_banded : eigenvalues and right eigenvectors for symmetric/Hermitian
  155. band matrices
  156. eigh_tridiagonal : eigenvalues and right eigenvectors for
  157. symmetric/Hermitian tridiagonal matrices
  158. Examples
  159. --------
  160. >>> import numpy as np
  161. >>> from scipy import linalg
  162. >>> a = np.array([[0., -1.], [1., 0.]])
  163. >>> linalg.eigvals(a)
  164. array([0.+1.j, 0.-1.j])
  165. >>> b = np.array([[0., 1.], [1., 1.]])
  166. >>> linalg.eigvals(a, b)
  167. array([ 1.+0.j, -1.+0.j])
  168. >>> a = np.array([[3., 0., 0.], [0., 8., 0.], [0., 0., 7.]])
  169. >>> linalg.eigvals(a, homogeneous_eigvals=True)
  170. array([[3.+0.j, 8.+0.j, 7.+0.j],
  171. [1.+0.j, 1.+0.j, 1.+0.j]])
  172. >>> a = np.array([[0., -1.], [1., 0.]])
  173. >>> linalg.eigvals(a) == linalg.eig(a)[0]
  174. array([ True, True])
  175. >>> linalg.eig(a, left=True, right=False)[1] # normalized left eigenvector
  176. array([[-0.70710678+0.j , -0.70710678-0.j ],
  177. [-0. +0.70710678j, -0. -0.70710678j]])
  178. >>> linalg.eig(a, left=False, right=True)[1] # normalized right eigenvector
  179. array([[0.70710678+0.j , 0.70710678-0.j ],
  180. [0. -0.70710678j, 0. +0.70710678j]])
  181. """
  182. a1 = _asarray_validated(a, check_finite=check_finite)
  183. if len(a1.shape) != 2 or a1.shape[0] != a1.shape[1]:
  184. raise ValueError('expected square matrix')
  185. # accommodate square empty matrices
  186. if a1.size == 0:
  187. w_n, vr_n = eig(np.eye(2, dtype=a1.dtype))
  188. w = np.empty_like(a1, shape=(0,), dtype=w_n.dtype)
  189. w = _make_eigvals(w, None, homogeneous_eigvals)
  190. vl = np.empty_like(a1, shape=(0, 0), dtype=vr_n.dtype)
  191. vr = np.empty_like(a1, shape=(0, 0), dtype=vr_n.dtype)
  192. if not (left or right):
  193. return w
  194. if left:
  195. if right:
  196. return w, vl, vr
  197. return w, vl
  198. return w, vr
  199. overwrite_a = overwrite_a or (_datacopied(a1, a))
  200. if b is not None:
  201. b1 = _asarray_validated(b, check_finite=check_finite)
  202. overwrite_b = overwrite_b or _datacopied(b1, b)
  203. if len(b1.shape) != 2 or b1.shape[0] != b1.shape[1]:
  204. raise ValueError('expected square matrix')
  205. if b1.shape != a1.shape:
  206. raise ValueError('a and b must have the same shape')
  207. return _geneig(a1, b1, left, right, overwrite_a, overwrite_b,
  208. homogeneous_eigvals)
  209. geev, geev_lwork = get_lapack_funcs(('geev', 'geev_lwork'), (a1,))
  210. compute_vl, compute_vr = left, right
  211. lwork = _compute_lwork(geev_lwork, a1.shape[0],
  212. compute_vl=compute_vl,
  213. compute_vr=compute_vr)
  214. if geev.typecode in 'cz':
  215. w, vl, vr, info = geev(a1, lwork=lwork,
  216. compute_vl=compute_vl,
  217. compute_vr=compute_vr,
  218. overwrite_a=overwrite_a)
  219. w = _make_eigvals(w, None, homogeneous_eigvals)
  220. else:
  221. wr, wi, vl, vr, info = geev(a1, lwork=lwork,
  222. compute_vl=compute_vl,
  223. compute_vr=compute_vr,
  224. overwrite_a=overwrite_a)
  225. w = wr + _I * wi
  226. w = _make_eigvals(w, None, homogeneous_eigvals)
  227. _check_info(info, 'eig algorithm (geev)',
  228. positive='did not converge (only eigenvalues '
  229. 'with order >= %d have converged)')
  230. only_real = np.all(w.imag == 0.0)
  231. if not (geev.typecode in 'cz' or only_real):
  232. t = w.dtype.char
  233. if left:
  234. vl = _make_complex_eigvecs(w, vl, t)
  235. if right:
  236. vr = _make_complex_eigvecs(w, vr, t)
  237. if not (left or right):
  238. return w
  239. if left:
  240. if right:
  241. return w, vl, vr
  242. return w, vl
  243. return w, vr
  244. @_apply_over_batch(('a', 2), ('b', 2))
  245. def eigh(a, b=None, *, lower=True, eigvals_only=False, overwrite_a=False,
  246. overwrite_b=False, type=1, check_finite=True, subset_by_index=None,
  247. subset_by_value=None, driver=None):
  248. """
  249. Solve a standard or generalized eigenvalue problem for a complex
  250. Hermitian or real symmetric matrix.
  251. Find eigenvalues array ``w`` and optionally eigenvectors array ``v`` of
  252. array ``a``, where ``b`` is positive definite such that for every
  253. eigenvalue λ (i-th entry of w) and its eigenvector ``vi`` (i-th column of
  254. ``v``) satisfies::
  255. a @ vi = λ * b @ vi
  256. vi.conj().T @ a @ vi = λ
  257. vi.conj().T @ b @ vi = 1
  258. In the standard problem, ``b`` is assumed to be the identity matrix.
  259. Parameters
  260. ----------
  261. a : (M, M) array_like
  262. A complex Hermitian or real symmetric matrix whose eigenvalues and
  263. eigenvectors will be computed.
  264. b : (M, M) array_like, optional
  265. A complex Hermitian or real symmetric definite positive matrix in.
  266. If omitted, identity matrix is assumed.
  267. lower : bool, optional
  268. Whether the pertinent array data is taken from the lower or upper
  269. triangle of ``a`` and, if applicable, ``b``. (Default: lower)
  270. eigvals_only : bool, optional
  271. Whether to calculate only eigenvalues and no eigenvectors.
  272. (Default: both are calculated)
  273. subset_by_index : iterable, optional
  274. If provided, this two-element iterable defines the start and the end
  275. indices of the desired eigenvalues (ascending order and 0-indexed).
  276. To return only the second smallest to fifth smallest eigenvalues,
  277. ``[1, 4]`` is used. ``[n-3, n-1]`` returns the largest three. Only
  278. available with "evr", "evx", and "gvx" drivers. The entries are
  279. directly converted to integers via ``int()``.
  280. subset_by_value : iterable, optional
  281. If provided, this two-element iterable defines the half-open interval
  282. ``(a, b]`` that, if any, only the eigenvalues between these values
  283. are returned. Only available with "evr", "evx", and "gvx" drivers. Use
  284. ``np.inf`` for the unconstrained ends.
  285. driver : str, optional
  286. Defines which LAPACK driver should be used. Valid options are "ev",
  287. "evd", "evr", "evx" for standard problems and "gv", "gvd", "gvx" for
  288. generalized (where b is not None) problems. See the Notes section.
  289. The default for standard problems is "evr". For generalized problems,
  290. "gvd" is used for full set, and "gvx" for subset requested cases.
  291. type : int, optional
  292. For the generalized problems, this keyword specifies the problem type
  293. to be solved for ``w`` and ``v`` (only takes 1, 2, 3 as possible
  294. inputs)::
  295. 1 => a @ v = w @ b @ v
  296. 2 => a @ b @ v = w @ v
  297. 3 => b @ a @ v = w @ v
  298. This keyword is ignored for standard problems.
  299. overwrite_a : bool, optional
  300. Whether to overwrite data in ``a`` (may improve performance). Default
  301. is False.
  302. overwrite_b : bool, optional
  303. Whether to overwrite data in ``b`` (may improve performance). Default
  304. is False.
  305. check_finite : bool, optional
  306. Whether to check that the input matrices contain only finite numbers.
  307. Disabling may give a performance gain, but may result in problems
  308. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  309. Returns
  310. -------
  311. w : (N,) ndarray
  312. The N (N<=M) selected eigenvalues, in ascending order, each
  313. repeated according to its multiplicity.
  314. v : (M, N) ndarray
  315. The normalized eigenvector corresponding to the eigenvalue ``w[i]`` is
  316. the column ``v[:,i]``. Only returned if ``eigvals_only=False``.
  317. Raises
  318. ------
  319. LinAlgError
  320. If eigenvalue computation does not converge, an error occurred, or
  321. b matrix is not definite positive. Note that if input matrices are
  322. not symmetric or Hermitian, no error will be reported but results will
  323. be wrong.
  324. See Also
  325. --------
  326. eigvalsh : eigenvalues of symmetric or Hermitian arrays
  327. eig : eigenvalues and right eigenvectors for non-symmetric arrays
  328. eigh_tridiagonal : eigenvalues and right eigenvectors for
  329. symmetric/Hermitian tridiagonal matrices
  330. Notes
  331. -----
  332. This function does not check the input array for being Hermitian/symmetric
  333. in order to allow for representing arrays with only their upper/lower
  334. triangular parts. Also, note that even though not taken into account,
  335. finiteness check applies to the whole array and unaffected by "lower"
  336. keyword.
  337. This function uses LAPACK drivers for computations in all possible keyword
  338. combinations, prefixed with ``sy`` if arrays are real and ``he`` if
  339. complex, e.g., a float array with "evr" driver is solved via
  340. "syevr", complex arrays with "gvx" driver problem is solved via "hegvx"
  341. etc.
  342. As a brief summary, the slowest and the most robust driver is the
  343. classical ``<sy/he>ev`` which uses symmetric QR. ``<sy/he>evr`` is seen as
  344. the optimal choice for the most general cases. However, there are certain
  345. occasions that ``<sy/he>evd`` computes faster at the expense of more
  346. memory usage. ``<sy/he>evx``, while still being faster than ``<sy/he>ev``,
  347. often performs worse than the rest except when very few eigenvalues are
  348. requested for large arrays though there is still no performance guarantee.
  349. Note that the underlying LAPACK algorithms are different depending on whether
  350. `eigvals_only` is True or False --- thus the eigenvalues may differ
  351. depending on whether eigenvectors are requested or not. The difference is
  352. generally of the order of machine epsilon times the largest eigenvalue,
  353. so is likely only visible for zero or nearly zero eigenvalues.
  354. For the generalized problem, normalization with respect to the given
  355. type argument::
  356. type 1 and 3 : v.conj().T @ a @ v = w
  357. type 2 : inv(v).conj().T @ a @ inv(v) = w
  358. type 1 or 2 : v.conj().T @ b @ v = I
  359. type 3 : v.conj().T @ inv(b) @ v = I
  360. Examples
  361. --------
  362. >>> import numpy as np
  363. >>> from scipy.linalg import eigh
  364. >>> A = np.array([[6, 3, 1, 5], [3, 0, 5, 1], [1, 5, 6, 2], [5, 1, 2, 2]])
  365. >>> w, v = eigh(A)
  366. >>> np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4)))
  367. True
  368. Request only the eigenvalues
  369. >>> w = eigh(A, eigvals_only=True)
  370. Request eigenvalues that are less than 10.
  371. >>> A = np.array([[34, -4, -10, -7, 2],
  372. ... [-4, 7, 2, 12, 0],
  373. ... [-10, 2, 44, 2, -19],
  374. ... [-7, 12, 2, 79, -34],
  375. ... [2, 0, -19, -34, 29]])
  376. >>> eigh(A, eigvals_only=True, subset_by_value=[-np.inf, 10])
  377. array([6.69199443e-07, 9.11938152e+00])
  378. Request the second smallest eigenvalue and its eigenvector
  379. >>> w, v = eigh(A, subset_by_index=[1, 1])
  380. >>> w
  381. array([9.11938152])
  382. >>> v.shape # only a single column is returned
  383. (5, 1)
  384. """
  385. # set lower
  386. uplo = 'L' if lower else 'U'
  387. # Set job for Fortran routines
  388. _job = 'N' if eigvals_only else 'V'
  389. drv_str = [None, "ev", "evd", "evr", "evx", "gv", "gvd", "gvx"]
  390. if driver not in drv_str:
  391. raise ValueError('"{}" is unknown. Possible values are "None", "{}".'
  392. ''.format(driver, '", "'.join(drv_str[1:])))
  393. a1 = _asarray_validated(a, check_finite=check_finite)
  394. if len(a1.shape) != 2 or a1.shape[0] != a1.shape[1]:
  395. raise ValueError('expected square "a" matrix')
  396. # accommodate square empty matrices
  397. if a1.size == 0:
  398. w_n, v_n = eigh(np.eye(2, dtype=a1.dtype))
  399. w = np.empty_like(a1, shape=(0,), dtype=w_n.dtype)
  400. v = np.empty_like(a1, shape=(0, 0), dtype=v_n.dtype)
  401. if eigvals_only:
  402. return w
  403. else:
  404. return w, v
  405. overwrite_a = overwrite_a or (_datacopied(a1, a))
  406. cplx = True if iscomplexobj(a1) else False
  407. n = a1.shape[0]
  408. drv_args = {'overwrite_a': overwrite_a}
  409. if b is not None:
  410. b1 = _asarray_validated(b, check_finite=check_finite)
  411. overwrite_b = overwrite_b or _datacopied(b1, b)
  412. if len(b1.shape) != 2 or b1.shape[0] != b1.shape[1]:
  413. raise ValueError('expected square "b" matrix')
  414. if b1.shape != a1.shape:
  415. raise ValueError(f"wrong b dimensions {b1.shape}, should be {a1.shape}")
  416. if type not in [1, 2, 3]:
  417. raise ValueError('"type" keyword only accepts 1, 2, and 3.')
  418. cplx = True if iscomplexobj(b1) else (cplx or False)
  419. drv_args.update({'overwrite_b': overwrite_b, 'itype': type})
  420. subset = (subset_by_index is not None) or (subset_by_value is not None)
  421. # Both subsets can't be given
  422. if subset_by_index and subset_by_value:
  423. raise ValueError('Either index or value subset can be requested.')
  424. # Check indices if given
  425. if subset_by_index:
  426. lo, hi = (int(x) for x in subset_by_index)
  427. if not (0 <= lo <= hi < n):
  428. raise ValueError('Requested eigenvalue indices are not valid. '
  429. f'Valid range is [0, {n-1}] and start <= end, but '
  430. f'start={lo}, end={hi} is given')
  431. # fortran is 1-indexed
  432. drv_args.update({'range': 'I', 'il': lo + 1, 'iu': hi + 1})
  433. if subset_by_value:
  434. lo, hi = subset_by_value
  435. if not (-inf <= lo < hi <= inf):
  436. raise ValueError('Requested eigenvalue bounds are not valid. '
  437. 'Valid range is (-inf, inf) and low < high, but '
  438. f'low={lo}, high={hi} is given')
  439. drv_args.update({'range': 'V', 'vl': lo, 'vu': hi})
  440. # fix prefix for lapack routines
  441. pfx = 'he' if cplx else 'sy'
  442. # decide on the driver if not given
  443. # first early exit on incompatible choice
  444. if driver:
  445. if b is None and (driver in ["gv", "gvd", "gvx"]):
  446. raise ValueError(f'{driver} requires input b array to be supplied '
  447. 'for generalized eigenvalue problems.')
  448. if (b is not None) and (driver in ['ev', 'evd', 'evr', 'evx']):
  449. raise ValueError(f'"{driver}" does not accept input b array '
  450. 'for standard eigenvalue problems.')
  451. if subset and (driver in ["ev", "evd", "gv", "gvd"]):
  452. raise ValueError(f'"{driver}" cannot compute subsets of eigenvalues')
  453. # Default driver is evr and gvd
  454. else:
  455. driver = "evr" if b is None else ("gvx" if subset else "gvd")
  456. lwork_spec = {
  457. 'syevd': ['lwork', 'liwork'],
  458. 'syevr': ['lwork', 'liwork'],
  459. 'heevd': ['lwork', 'liwork', 'lrwork'],
  460. 'heevr': ['lwork', 'lrwork', 'liwork'],
  461. }
  462. if b is None: # Standard problem
  463. drv, drvlw = get_lapack_funcs((pfx + driver, pfx+driver+'_lwork'),
  464. [a1])
  465. clw_args = {'n': n, 'lower': lower}
  466. if driver == 'evd':
  467. clw_args.update({'compute_v': 0 if _job == "N" else 1})
  468. lw = _compute_lwork(drvlw, **clw_args)
  469. # Multiple lwork vars
  470. if isinstance(lw, tuple):
  471. lwork_args = dict(zip(lwork_spec[pfx+driver], lw))
  472. else:
  473. lwork_args = {'lwork': lw}
  474. drv_args.update({'lower': lower, 'compute_v': 0 if _job == "N" else 1})
  475. w, v, *other_args, info = drv(a=a1, **drv_args, **lwork_args)
  476. else: # Generalized problem
  477. # 'gvd' doesn't have lwork query
  478. if driver == "gvd":
  479. drv = get_lapack_funcs(pfx + "gvd", [a1, b1])
  480. lwork_args = {}
  481. else:
  482. drv, drvlw = get_lapack_funcs((pfx + driver, pfx+driver+'_lwork'),
  483. [a1, b1])
  484. # generalized drivers use uplo instead of lower
  485. lw = _compute_lwork(drvlw, n, uplo=uplo)
  486. lwork_args = {'lwork': lw}
  487. drv_args.update({'uplo': uplo, 'jobz': _job})
  488. w, v, *other_args, info = drv(a=a1, b=b1, **drv_args, **lwork_args)
  489. # m is always the first extra argument
  490. w = w[:other_args[0]] if subset else w
  491. v = v[:, :other_args[0]] if (subset and not eigvals_only) else v
  492. # Check if we had a successful exit
  493. if info == 0:
  494. if eigvals_only:
  495. return w
  496. else:
  497. return w, v
  498. else:
  499. if info < -1:
  500. raise LinAlgError(f'Illegal value in argument {-info} of internal '
  501. f'{drv.typecode + pfx + driver}')
  502. elif info > n:
  503. raise LinAlgError(f'The leading minor of order {info-n} of B is not '
  504. 'positive definite. The factorization of B '
  505. 'could not be completed and no eigenvalues '
  506. 'or eigenvectors were computed.')
  507. else:
  508. drv_err = {'ev': 'The algorithm failed to converge; {} '
  509. 'off-diagonal elements of an intermediate '
  510. 'tridiagonal form did not converge to zero.',
  511. 'evx': '{} eigenvectors failed to converge.',
  512. 'evd': 'The algorithm failed to compute an eigenvalue '
  513. 'while working on the submatrix lying in rows '
  514. 'and columns {0}/{1} through mod({0},{1}).',
  515. 'evr': 'Internal Error.'
  516. }
  517. if driver in ['ev', 'gv']:
  518. msg = drv_err['ev'].format(info)
  519. elif driver in ['evx', 'gvx']:
  520. msg = drv_err['evx'].format(info)
  521. elif driver in ['evd', 'gvd']:
  522. if eigvals_only:
  523. msg = drv_err['ev'].format(info)
  524. else:
  525. msg = drv_err['evd'].format(info, n+1)
  526. else:
  527. msg = drv_err['evr']
  528. raise LinAlgError(msg)
  529. _conv_dict = {0: 0, 1: 1, 2: 2,
  530. 'all': 0, 'value': 1, 'index': 2,
  531. 'a': 0, 'v': 1, 'i': 2}
  532. def _check_select(select, select_range, max_ev, max_len):
  533. """Check that select is valid, convert to Fortran style."""
  534. if isinstance(select, str):
  535. select = select.lower()
  536. try:
  537. select = _conv_dict[select]
  538. except KeyError as e:
  539. raise ValueError('invalid argument for select') from e
  540. vl, vu = 0., 1.
  541. il = iu = 1
  542. if select != 0: # (non-all)
  543. sr = asarray(select_range)
  544. if sr.ndim != 1 or sr.size != 2 or sr[1] < sr[0]:
  545. raise ValueError('select_range must be a 2-element array-like '
  546. 'in nondecreasing order')
  547. if select == 1: # (value)
  548. vl, vu = sr
  549. if max_ev == 0:
  550. max_ev = max_len
  551. else: # 2 (index)
  552. if sr.dtype.char.lower() not in 'hilqp':
  553. raise ValueError(
  554. f'when using select="i", select_range must '
  555. f'contain integers, got dtype {sr.dtype} ({sr.dtype.char})'
  556. )
  557. # translate Python (0 ... N-1) into Fortran (1 ... N) with + 1
  558. il, iu = sr + 1
  559. if min(il, iu) < 1 or max(il, iu) > max_len:
  560. raise ValueError('select_range out of bounds')
  561. max_ev = iu - il + 1
  562. return select, vl, vu, il, iu, max_ev
  563. @_apply_over_batch(('a_band', 2))
  564. def eig_banded(a_band, lower=False, eigvals_only=False, overwrite_a_band=False,
  565. select='a', select_range=None, max_ev=0, check_finite=True):
  566. """
  567. Solve real symmetric or complex Hermitian band matrix eigenvalue problem.
  568. Find eigenvalues w and optionally right eigenvectors v of a::
  569. a v[:,i] = w[i] v[:,i]
  570. v.H v = identity
  571. The matrix a is stored in a_band either in lower diagonal or upper
  572. diagonal ordered form:
  573. a_band[u + i - j, j] == a[i,j] (if upper form; i <= j)
  574. a_band[ i - j, j] == a[i,j] (if lower form; i >= j)
  575. where u is the number of bands above the diagonal.
  576. Example of a_band (shape of a is (6,6), u=2)::
  577. upper form:
  578. * * a02 a13 a24 a35
  579. * a01 a12 a23 a34 a45
  580. a00 a11 a22 a33 a44 a55
  581. lower form:
  582. a00 a11 a22 a33 a44 a55
  583. a10 a21 a32 a43 a54 *
  584. a20 a31 a42 a53 * *
  585. Cells marked with * are not used.
  586. Parameters
  587. ----------
  588. a_band : (u+1, M) array_like
  589. The bands of the M by M matrix a.
  590. lower : bool, optional
  591. Is the matrix in the lower form. (Default is upper form)
  592. eigvals_only : bool, optional
  593. Compute only the eigenvalues and no eigenvectors.
  594. (Default: calculate also eigenvectors)
  595. overwrite_a_band : bool, optional
  596. Discard data in a_band (may enhance performance)
  597. select : {'a', 'v', 'i'}, optional
  598. Which eigenvalues to calculate
  599. ====== ========================================
  600. select calculated
  601. ====== ========================================
  602. 'a' All eigenvalues
  603. 'v' Eigenvalues in the interval (min, max]
  604. 'i' Eigenvalues with indices min <= i <= max
  605. ====== ========================================
  606. select_range : (min, max), optional
  607. Range of selected eigenvalues
  608. max_ev : int, optional
  609. For select=='v', maximum number of eigenvalues expected.
  610. For other values of select, has no meaning.
  611. In doubt, leave this parameter untouched.
  612. check_finite : bool, optional
  613. Whether to check that the input matrix contains only finite numbers.
  614. Disabling may give a performance gain, but may result in problems
  615. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  616. Returns
  617. -------
  618. w : (M,) ndarray
  619. The eigenvalues, in ascending order, each repeated according to its
  620. multiplicity.
  621. v : (M, M) float or complex ndarray
  622. The normalized eigenvector corresponding to the eigenvalue w[i] is
  623. the column v[:,i]. Only returned if ``eigvals_only=False``.
  624. Raises
  625. ------
  626. LinAlgError
  627. If eigenvalue computation does not converge.
  628. See Also
  629. --------
  630. eigvals_banded : eigenvalues for symmetric/Hermitian band matrices
  631. eig : eigenvalues and right eigenvectors of general arrays.
  632. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  633. eigh_tridiagonal : eigenvalues and right eigenvectors for
  634. symmetric/Hermitian tridiagonal matrices
  635. Examples
  636. --------
  637. >>> import numpy as np
  638. >>> from scipy.linalg import eig_banded
  639. >>> A = np.array([[1, 5, 2, 0], [5, 2, 5, 2], [2, 5, 3, 5], [0, 2, 5, 4]])
  640. >>> Ab = np.array([[1, 2, 3, 4], [5, 5, 5, 0], [2, 2, 0, 0]])
  641. >>> w, v = eig_banded(Ab, lower=True)
  642. >>> np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4)))
  643. True
  644. >>> w = eig_banded(Ab, lower=True, eigvals_only=True)
  645. >>> w
  646. array([-4.26200532, -2.22987175, 3.95222349, 12.53965359])
  647. Request only the eigenvalues between ``[-3, 4]``
  648. >>> w, v = eig_banded(Ab, lower=True, select='v', select_range=[-3, 4])
  649. >>> w
  650. array([-2.22987175, 3.95222349])
  651. """
  652. if eigvals_only or overwrite_a_band:
  653. a1 = _asarray_validated(a_band, check_finite=check_finite)
  654. overwrite_a_band = overwrite_a_band or (_datacopied(a1, a_band))
  655. else:
  656. a1 = array(a_band)
  657. if issubclass(a1.dtype.type, inexact) and not isfinite(a1).all():
  658. raise ValueError("array must not contain infs or NaNs")
  659. overwrite_a_band = 1
  660. if len(a1.shape) != 2:
  661. raise ValueError('expected a 2-D array')
  662. # accommodate square empty matrices
  663. if a1.size == 0:
  664. w_n, v_n = eig_banded(np.array([[0, 0], [1, 1]], dtype=a1.dtype))
  665. w = np.empty_like(a1, shape=(0,), dtype=w_n.dtype)
  666. v = np.empty_like(a1, shape=(0, 0), dtype=v_n.dtype)
  667. if eigvals_only:
  668. return w
  669. else:
  670. return w, v
  671. select, vl, vu, il, iu, max_ev = _check_select(
  672. select, select_range, max_ev, a1.shape[1])
  673. del select_range
  674. if select == 0:
  675. if a1.dtype.char in 'GFD':
  676. # FIXME: implement this somewhen, for now go with builtin values
  677. # FIXME: calc optimal lwork by calling ?hbevd(lwork=-1)
  678. # or by using calc_lwork.f ???
  679. # lwork = calc_lwork.hbevd(bevd.typecode, a1.shape[0], lower)
  680. internal_name = 'hbevd'
  681. else: # a1.dtype.char in 'fd':
  682. # FIXME: implement this somewhen, for now go with builtin values
  683. # see above
  684. # lwork = calc_lwork.sbevd(bevd.typecode, a1.shape[0], lower)
  685. internal_name = 'sbevd'
  686. bevd, = get_lapack_funcs((internal_name,), (a1,))
  687. w, v, info = bevd(a1, compute_v=not eigvals_only,
  688. lower=lower, overwrite_ab=overwrite_a_band)
  689. else: # select in [1, 2]
  690. if eigvals_only:
  691. max_ev = 1
  692. # calculate optimal abstol for dsbevx (see manpage)
  693. if a1.dtype.char in 'fF': # single precision
  694. lamch, = get_lapack_funcs(('lamch',), (array(0, dtype='f'),))
  695. else:
  696. lamch, = get_lapack_funcs(('lamch',), (array(0, dtype='d'),))
  697. abstol = 2 * lamch('s')
  698. if a1.dtype.char in 'GFD':
  699. internal_name = 'hbevx'
  700. else: # a1.dtype.char in 'gfd'
  701. internal_name = 'sbevx'
  702. bevx, = get_lapack_funcs((internal_name,), (a1,))
  703. w, v, m, ifail, info = bevx(
  704. a1, vl, vu, il, iu, compute_v=not eigvals_only, mmax=max_ev,
  705. range=select, lower=lower, overwrite_ab=overwrite_a_band,
  706. abstol=abstol)
  707. # crop off w and v
  708. w = w[:m]
  709. if not eigvals_only:
  710. v = v[:, :m]
  711. _check_info(info, internal_name)
  712. if eigvals_only:
  713. return w
  714. return w, v
  715. @_apply_over_batch(('a', 2), ('b', 2))
  716. def eigvals(a, b=None, overwrite_a=False, check_finite=True,
  717. homogeneous_eigvals=False):
  718. """
  719. Compute eigenvalues from an ordinary or generalized eigenvalue problem.
  720. Find eigenvalues of a general matrix::
  721. a vr[:,i] = w[i] b vr[:,i]
  722. Parameters
  723. ----------
  724. a : (M, M) array_like
  725. A complex or real matrix whose eigenvalues and eigenvectors
  726. will be computed.
  727. b : (M, M) array_like, optional
  728. Right-hand side matrix in a generalized eigenvalue problem.
  729. If omitted, identity matrix is assumed.
  730. overwrite_a : bool, optional
  731. Whether to overwrite data in a (may improve performance)
  732. check_finite : bool, optional
  733. Whether to check that the input matrices contain only finite numbers.
  734. Disabling may give a performance gain, but may result in problems
  735. (crashes, non-termination) if the inputs do contain infinities
  736. or NaNs.
  737. homogeneous_eigvals : bool, optional
  738. If True, return the eigenvalues in homogeneous coordinates.
  739. In this case ``w`` is a (2, M) array so that::
  740. w[1,i] a vr[:,i] = w[0,i] b vr[:,i]
  741. Default is False.
  742. Returns
  743. -------
  744. w : (M,) or (2, M) double or complex ndarray
  745. The eigenvalues, each repeated according to its multiplicity
  746. but not in any specific order. The shape is (M,) unless
  747. ``homogeneous_eigvals=True``.
  748. Raises
  749. ------
  750. LinAlgError
  751. If eigenvalue computation does not converge
  752. See Also
  753. --------
  754. eig : eigenvalues and right eigenvectors of general arrays.
  755. eigvalsh : eigenvalues of symmetric or Hermitian arrays
  756. eigvals_banded : eigenvalues for symmetric/Hermitian band matrices
  757. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  758. matrices
  759. Examples
  760. --------
  761. >>> import numpy as np
  762. >>> from scipy import linalg
  763. >>> a = np.array([[0., -1.], [1., 0.]])
  764. >>> linalg.eigvals(a)
  765. array([0.+1.j, 0.-1.j])
  766. >>> b = np.array([[0., 1.], [1., 1.]])
  767. >>> linalg.eigvals(a, b)
  768. array([ 1.+0.j, -1.+0.j])
  769. >>> a = np.array([[3., 0., 0.], [0., 8., 0.], [0., 0., 7.]])
  770. >>> linalg.eigvals(a, homogeneous_eigvals=True)
  771. array([[3.+0.j, 8.+0.j, 7.+0.j],
  772. [1.+0.j, 1.+0.j, 1.+0.j]])
  773. """
  774. return eig(a, b=b, left=0, right=0, overwrite_a=overwrite_a,
  775. check_finite=check_finite,
  776. homogeneous_eigvals=homogeneous_eigvals)
  777. @_apply_over_batch(('a', 2), ('b', 2))
  778. def eigvalsh(a, b=None, *, lower=True, overwrite_a=False,
  779. overwrite_b=False, type=1, check_finite=True, subset_by_index=None,
  780. subset_by_value=None, driver=None):
  781. """
  782. Solves a standard or generalized eigenvalue problem for a complex
  783. Hermitian or real symmetric matrix.
  784. Find eigenvalues array ``w`` of array ``a``, where ``b`` is positive
  785. definite such that for every eigenvalue λ (i-th entry of w) and its
  786. eigenvector vi (i-th column of v) satisfies::
  787. a @ vi = λ * b @ vi
  788. vi.conj().T @ a @ vi = λ
  789. vi.conj().T @ b @ vi = 1
  790. In the standard problem, b is assumed to be the identity matrix.
  791. Parameters
  792. ----------
  793. a : (M, M) array_like
  794. A complex Hermitian or real symmetric matrix whose eigenvalues will
  795. be computed.
  796. b : (M, M) array_like, optional
  797. A complex Hermitian or real symmetric definite positive matrix in.
  798. If omitted, identity matrix is assumed.
  799. lower : bool, optional
  800. Whether the pertinent array data is taken from the lower or upper
  801. triangle of ``a`` and, if applicable, ``b``. (Default: lower)
  802. overwrite_a : bool, optional
  803. Whether to overwrite data in ``a`` (may improve performance). Default
  804. is False.
  805. overwrite_b : bool, optional
  806. Whether to overwrite data in ``b`` (may improve performance). Default
  807. is False.
  808. type : int, optional
  809. For the generalized problems, this keyword specifies the problem type
  810. to be solved for ``w`` and ``v`` (only takes 1, 2, 3 as possible
  811. inputs)::
  812. 1 => a @ v = w @ b @ v
  813. 2 => a @ b @ v = w @ v
  814. 3 => b @ a @ v = w @ v
  815. This keyword is ignored for standard problems.
  816. check_finite : bool, optional
  817. Whether to check that the input matrices contain only finite numbers.
  818. Disabling may give a performance gain, but may result in problems
  819. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  820. subset_by_index : iterable, optional
  821. If provided, this two-element iterable defines the start and the end
  822. indices of the desired eigenvalues (ascending order and 0-indexed).
  823. To return only the second smallest to fifth smallest eigenvalues,
  824. ``[1, 4]`` is used. ``[n-3, n-1]`` returns the largest three. Only
  825. available with "evr", "evx", and "gvx" drivers. The entries are
  826. directly converted to integers via ``int()``.
  827. subset_by_value : iterable, optional
  828. If provided, this two-element iterable defines the half-open interval
  829. ``(a, b]`` that, if any, only the eigenvalues between these values
  830. are returned. Only available with "evr", "evx", and "gvx" drivers. Use
  831. ``np.inf`` for the unconstrained ends.
  832. driver : str, optional
  833. Defines which LAPACK driver should be used. Valid options are "ev",
  834. "evd", "evr", "evx" for standard problems and "gv", "gvd", "gvx" for
  835. generalized (where b is not None) problems. See the Notes section of
  836. `scipy.linalg.eigh`.
  837. Returns
  838. -------
  839. w : (N,) ndarray
  840. The N (N<=M) selected eigenvalues, in ascending order, each
  841. repeated according to its multiplicity.
  842. Raises
  843. ------
  844. LinAlgError
  845. If eigenvalue computation does not converge, an error occurred, or
  846. b matrix is not definite positive. Note that if input matrices are
  847. not symmetric or Hermitian, no error will be reported but results will
  848. be wrong.
  849. See Also
  850. --------
  851. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  852. eigvals : eigenvalues of general arrays
  853. eigvals_banded : eigenvalues for symmetric/Hermitian band matrices
  854. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  855. matrices
  856. Notes
  857. -----
  858. This function does not check the input array for being Hermitian/symmetric
  859. in order to allow for representing arrays with only their upper/lower
  860. triangular parts.
  861. This function serves as a one-liner shorthand for `scipy.linalg.eigh` with
  862. the option ``eigvals_only=True`` to get the eigenvalues and not the
  863. eigenvectors. Here it is kept as a legacy convenience. It might be
  864. beneficial to use the main function to have full control and to be a bit
  865. more pythonic.
  866. Examples
  867. --------
  868. For more examples see `scipy.linalg.eigh`.
  869. >>> import numpy as np
  870. >>> from scipy.linalg import eigvalsh
  871. >>> A = np.array([[6, 3, 1, 5], [3, 0, 5, 1], [1, 5, 6, 2], [5, 1, 2, 2]])
  872. >>> w = eigvalsh(A)
  873. >>> w
  874. array([-3.74637491, -0.76263923, 6.08502336, 12.42399079])
  875. """
  876. return eigh(a, b=b, lower=lower, eigvals_only=True, overwrite_a=overwrite_a,
  877. overwrite_b=overwrite_b, type=type, check_finite=check_finite,
  878. subset_by_index=subset_by_index, subset_by_value=subset_by_value,
  879. driver=driver)
  880. @_apply_over_batch(('a_band', 2))
  881. def eigvals_banded(a_band, lower=False, overwrite_a_band=False,
  882. select='a', select_range=None, check_finite=True):
  883. """
  884. Solve real symmetric or complex Hermitian band matrix eigenvalue problem.
  885. Find eigenvalues w of a::
  886. a v[:,i] = w[i] v[:,i]
  887. v.H v = identity
  888. The matrix a is stored in a_band either in lower diagonal or upper
  889. diagonal ordered form:
  890. a_band[u + i - j, j] == a[i,j] (if upper form; i <= j)
  891. a_band[ i - j, j] == a[i,j] (if lower form; i >= j)
  892. where u is the number of bands above the diagonal.
  893. Example of a_band (shape of a is (6,6), u=2)::
  894. upper form:
  895. * * a02 a13 a24 a35
  896. * a01 a12 a23 a34 a45
  897. a00 a11 a22 a33 a44 a55
  898. lower form:
  899. a00 a11 a22 a33 a44 a55
  900. a10 a21 a32 a43 a54 *
  901. a20 a31 a42 a53 * *
  902. Cells marked with * are not used.
  903. Parameters
  904. ----------
  905. a_band : (u+1, M) array_like
  906. The bands of the M by M matrix a.
  907. lower : bool, optional
  908. Is the matrix in the lower form. (Default is upper form)
  909. overwrite_a_band : bool, optional
  910. Discard data in a_band (may enhance performance)
  911. select : {'a', 'v', 'i'}, optional
  912. Which eigenvalues to calculate
  913. ====== ========================================
  914. select calculated
  915. ====== ========================================
  916. 'a' All eigenvalues
  917. 'v' Eigenvalues in the interval (min, max]
  918. 'i' Eigenvalues with indices min <= i <= max
  919. ====== ========================================
  920. select_range : (min, max), optional
  921. Range of selected eigenvalues
  922. check_finite : bool, optional
  923. Whether to check that the input matrix contains only finite numbers.
  924. Disabling may give a performance gain, but may result in problems
  925. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  926. Returns
  927. -------
  928. w : (M,) ndarray
  929. The eigenvalues, in ascending order, each repeated according to its
  930. multiplicity.
  931. Raises
  932. ------
  933. LinAlgError
  934. If eigenvalue computation does not converge.
  935. See Also
  936. --------
  937. eig_banded : eigenvalues and right eigenvectors for symmetric/Hermitian
  938. band matrices
  939. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  940. matrices
  941. eigvals : eigenvalues of general arrays
  942. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  943. eig : eigenvalues and right eigenvectors for non-symmetric arrays
  944. Examples
  945. --------
  946. >>> import numpy as np
  947. >>> from scipy.linalg import eigvals_banded
  948. >>> A = np.array([[1, 5, 2, 0], [5, 2, 5, 2], [2, 5, 3, 5], [0, 2, 5, 4]])
  949. >>> Ab = np.array([[1, 2, 3, 4], [5, 5, 5, 0], [2, 2, 0, 0]])
  950. >>> w = eigvals_banded(Ab, lower=True)
  951. >>> w
  952. array([-4.26200532, -2.22987175, 3.95222349, 12.53965359])
  953. """
  954. return eig_banded(a_band, lower=lower, eigvals_only=1,
  955. overwrite_a_band=overwrite_a_band, select=select,
  956. select_range=select_range, check_finite=check_finite)
  957. @_apply_over_batch(('d', 1), ('e', 1))
  958. def eigvalsh_tridiagonal(d, e, select='a', select_range=None,
  959. check_finite=True, tol=0., lapack_driver='auto'):
  960. """
  961. Solve eigenvalue problem for a real symmetric tridiagonal matrix.
  962. Find eigenvalues `w` of ``a``::
  963. a v[:,i] = w[i] v[:,i]
  964. v.H v = identity
  965. For a real symmetric matrix ``a`` with diagonal elements `d` and
  966. off-diagonal elements `e`.
  967. Parameters
  968. ----------
  969. d : ndarray, shape (ndim,)
  970. The diagonal elements of the array.
  971. e : ndarray, shape (ndim-1,)
  972. The off-diagonal elements of the array.
  973. select : {'a', 'v', 'i'}, optional
  974. Which eigenvalues to calculate
  975. ====== ========================================
  976. select calculated
  977. ====== ========================================
  978. 'a' All eigenvalues
  979. 'v' Eigenvalues in the interval (min, max]
  980. 'i' Eigenvalues with indices min <= i <= max
  981. ====== ========================================
  982. select_range : (min, max), optional
  983. Range of selected eigenvalues
  984. check_finite : bool, optional
  985. Whether to check that the input matrix contains only finite numbers.
  986. Disabling may give a performance gain, but may result in problems
  987. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  988. tol : float
  989. The absolute tolerance to which each eigenvalue is required
  990. (only used when ``lapack_driver='stebz'``).
  991. An eigenvalue (or cluster) is considered to have converged if it
  992. lies in an interval of this width. If <= 0. (default),
  993. the value ``eps*|a|`` is used where eps is the machine precision,
  994. and ``|a|`` is the 1-norm of the matrix ``a``.
  995. lapack_driver : str
  996. LAPACK function to use, can be 'auto', 'stemr', 'stebz', 'sterf',
  997. 'stev', or 'stevd'. When 'auto' (default), it will use 'stevd' if
  998. ``select='a'`` and 'stebz' otherwise. 'sterf' and 'stev' can only
  999. be used when ``select='a'``.
  1000. Returns
  1001. -------
  1002. w : (M,) ndarray
  1003. The eigenvalues, in ascending order, each repeated according to its
  1004. multiplicity.
  1005. Raises
  1006. ------
  1007. LinAlgError
  1008. If eigenvalue computation does not converge.
  1009. See Also
  1010. --------
  1011. eigh_tridiagonal : eigenvalues and right eigenvectors for
  1012. symmetric/Hermitian tridiagonal matrices
  1013. Examples
  1014. --------
  1015. >>> import numpy as np
  1016. >>> from scipy.linalg import eigvalsh_tridiagonal, eigvalsh
  1017. >>> d = 3*np.ones(4)
  1018. >>> e = -1*np.ones(3)
  1019. >>> w = eigvalsh_tridiagonal(d, e)
  1020. >>> A = np.diag(d) + np.diag(e, k=1) + np.diag(e, k=-1)
  1021. >>> w2 = eigvalsh(A) # Verify with other eigenvalue routines
  1022. >>> np.allclose(w - w2, np.zeros(4))
  1023. True
  1024. """
  1025. return eigh_tridiagonal(
  1026. d, e, eigvals_only=True, select=select, select_range=select_range,
  1027. check_finite=check_finite, tol=tol, lapack_driver=lapack_driver)
  1028. @_apply_over_batch(('d', 1), ('e', 1))
  1029. def eigh_tridiagonal(d, e, eigvals_only=False, select='a', select_range=None,
  1030. check_finite=True, tol=0., lapack_driver='auto'):
  1031. """
  1032. Solve eigenvalue problem for a real symmetric tridiagonal matrix.
  1033. Find eigenvalues `w` and optionally right eigenvectors `v` of ``a``::
  1034. a v[:,i] = w[i] v[:,i]
  1035. v.H v = identity
  1036. For a real symmetric matrix ``a`` with diagonal elements `d` and
  1037. off-diagonal elements `e`.
  1038. Parameters
  1039. ----------
  1040. d : ndarray, shape (ndim,)
  1041. The diagonal elements of the array.
  1042. e : ndarray, shape (ndim-1,)
  1043. The off-diagonal elements of the array.
  1044. eigvals_only : bool, optional
  1045. Compute only the eigenvalues and no eigenvectors.
  1046. (Default: calculate also eigenvectors)
  1047. select : {'a', 'v', 'i'}, optional
  1048. Which eigenvalues to calculate
  1049. ====== ========================================
  1050. select calculated
  1051. ====== ========================================
  1052. 'a' All eigenvalues
  1053. 'v' Eigenvalues in the interval (min, max]
  1054. 'i' Eigenvalues with indices min <= i <= max
  1055. ====== ========================================
  1056. select_range : (min, max), optional
  1057. Range of selected eigenvalues
  1058. check_finite : bool, optional
  1059. Whether to check that the input matrix contains only finite numbers.
  1060. Disabling may give a performance gain, but may result in problems
  1061. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  1062. tol : float
  1063. The absolute tolerance to which each eigenvalue is required
  1064. (only used when 'stebz' is the `lapack_driver`).
  1065. An eigenvalue (or cluster) is considered to have converged if it
  1066. lies in an interval of this width. If <= 0. (default),
  1067. the value ``eps*|a|`` is used where eps is the machine precision,
  1068. and ``|a|`` is the 1-norm of the matrix ``a``.
  1069. lapack_driver : str
  1070. LAPACK function to use, can be 'auto', 'stemr', 'stebz', 'sterf',
  1071. 'stev', or 'stevd'. When 'auto' (default), it will use 'stevd' if ``select='a'``
  1072. and 'stebz' otherwise. When 'stebz' is used to find the eigenvalues and
  1073. ``eigvals_only=False``, then a second LAPACK call (to ``?STEIN``) is
  1074. used to find the corresponding eigenvectors. 'sterf' can only be
  1075. used when ``eigvals_only=True`` and ``select='a'``. 'stev' can only
  1076. be used when ``select='a'``.
  1077. Returns
  1078. -------
  1079. w : (M,) ndarray
  1080. The eigenvalues, in ascending order, each repeated according to its
  1081. multiplicity.
  1082. v : (M, M) ndarray
  1083. The normalized eigenvector corresponding to the eigenvalue ``w[i]`` is
  1084. the column ``v[:,i]``. Only returned if ``eigvals_only=False``.
  1085. Raises
  1086. ------
  1087. LinAlgError
  1088. If eigenvalue computation does not converge.
  1089. See Also
  1090. --------
  1091. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  1092. matrices
  1093. eig : eigenvalues and right eigenvectors for non-symmetric arrays
  1094. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  1095. eig_banded : eigenvalues and right eigenvectors for symmetric/Hermitian
  1096. band matrices
  1097. Notes
  1098. -----
  1099. This function makes use of LAPACK ``S/DSTEMR`` routines.
  1100. Examples
  1101. --------
  1102. >>> import numpy as np
  1103. >>> from scipy.linalg import eigh_tridiagonal
  1104. >>> d = 3*np.ones(4)
  1105. >>> e = -1*np.ones(3)
  1106. >>> w, v = eigh_tridiagonal(d, e)
  1107. >>> A = np.diag(d) + np.diag(e, k=1) + np.diag(e, k=-1)
  1108. >>> np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4)))
  1109. True
  1110. """
  1111. d = _asarray_validated(d, check_finite=check_finite)
  1112. e = _asarray_validated(e, check_finite=check_finite)
  1113. for check in (d, e):
  1114. if check.ndim != 1:
  1115. raise ValueError('expected a 1-D array')
  1116. if check.dtype.char in 'GFD': # complex
  1117. raise TypeError('Only real arrays currently supported')
  1118. if d.size != e.size + 1:
  1119. raise ValueError(f'd ({d.size}) must have one more element than e ({e.size})')
  1120. select, vl, vu, il, iu, _ = _check_select(
  1121. select, select_range, 0, d.size)
  1122. if not isinstance(lapack_driver, str):
  1123. raise TypeError('lapack_driver must be str')
  1124. drivers = ('auto', 'stemr', 'sterf', 'stebz', 'stev', 'stevd')
  1125. if lapack_driver not in drivers:
  1126. raise ValueError(f'lapack_driver must be one of {drivers}, '
  1127. f'got {lapack_driver}')
  1128. if lapack_driver == 'auto':
  1129. lapack_driver = 'stevd' if select == 0 else 'stebz'
  1130. # Quick exit for 1x1 case
  1131. if len(d) == 1:
  1132. if select == 1 and (not (vl < d[0] <= vu)): # request by value
  1133. w = array([])
  1134. v = empty([1, 0], dtype=d.dtype)
  1135. else: # all and request by index
  1136. w = array([d[0]], dtype=d.dtype)
  1137. v = array([[1.]], dtype=d.dtype)
  1138. if eigvals_only:
  1139. return w
  1140. else:
  1141. return w, v
  1142. func, = get_lapack_funcs((lapack_driver,), (d, e))
  1143. compute_v = not eigvals_only
  1144. if lapack_driver == 'sterf':
  1145. if select != 0:
  1146. raise ValueError('sterf can only be used when select == "a"')
  1147. if not eigvals_only:
  1148. raise ValueError('sterf can only be used when eigvals_only is '
  1149. 'True')
  1150. w, info = func(d, e)
  1151. m = len(w)
  1152. elif lapack_driver == 'stev':
  1153. if select != 0:
  1154. raise ValueError('stev can only be used when select == "a"')
  1155. w, v, info = func(d, e, compute_v=compute_v)
  1156. m = len(w)
  1157. elif lapack_driver == 'stevd':
  1158. if select != 0:
  1159. raise ValueError('stevd can only be used when select == "a"')
  1160. w, v, info = func(d, e, compute_v=compute_v)
  1161. m = len(w)
  1162. elif lapack_driver == 'stebz':
  1163. tol = float(tol)
  1164. internal_name = 'stebz'
  1165. stebz, = get_lapack_funcs((internal_name,), (d, e))
  1166. # If getting eigenvectors, needs to be block-ordered (B) instead of
  1167. # matrix-ordered (E), and we will reorder later
  1168. order = 'E' if eigvals_only else 'B'
  1169. m, w, iblock, isplit, info = stebz(d, e, select, vl, vu, il, iu, tol,
  1170. order)
  1171. else: # 'stemr'
  1172. # ?STEMR annoyingly requires size N instead of N-1
  1173. e_ = empty(e.size+1, e.dtype)
  1174. e_[:-1] = e
  1175. stemr_lwork, = get_lapack_funcs(('stemr_lwork',), (d, e))
  1176. lwork, liwork, info = stemr_lwork(d, e_, select, vl, vu, il, iu,
  1177. compute_v=compute_v)
  1178. _check_info(info, 'stemr_lwork')
  1179. m, w, v, info = func(d, e_, select, vl, vu, il, iu,
  1180. compute_v=compute_v, lwork=lwork, liwork=liwork)
  1181. _check_info(info, lapack_driver + ' (eigh_tridiagonal)')
  1182. w = w[:m]
  1183. if eigvals_only:
  1184. return w
  1185. else:
  1186. # Do we still need to compute the eigenvalues?
  1187. if lapack_driver == 'stebz':
  1188. func, = get_lapack_funcs(('stein',), (d, e))
  1189. v, info = func(d, e, w, iblock, isplit)
  1190. _check_info(info, 'stein (eigh_tridiagonal)',
  1191. positive='%d eigenvectors failed to converge')
  1192. # Convert block-order to matrix-order
  1193. order = argsort(w)
  1194. w, v = w[order], v[:, order]
  1195. else:
  1196. v = v[:, :m]
  1197. return w, v
  1198. def _check_info(info, driver, positive='did not converge (LAPACK info=%d)'):
  1199. """Check info return value."""
  1200. if info < 0:
  1201. raise ValueError(f'illegal value in argument {-info} of internal {driver}')
  1202. if info > 0 and positive:
  1203. raise LinAlgError(("%s " + positive) % (driver, info,))
  1204. @_apply_over_batch(('a', 2))
  1205. def hessenberg(a, calc_q=False, overwrite_a=False, check_finite=True):
  1206. """
  1207. Compute Hessenberg form of a matrix.
  1208. The Hessenberg decomposition is::
  1209. A = Q H Q^H
  1210. where `Q` is unitary/orthogonal and `H` has only zero elements below
  1211. the first sub-diagonal.
  1212. Parameters
  1213. ----------
  1214. a : (M, M) array_like
  1215. Matrix to bring into Hessenberg form.
  1216. calc_q : bool, optional
  1217. Whether to compute the transformation matrix. Default is False.
  1218. overwrite_a : bool, optional
  1219. Whether to overwrite `a`; may improve performance.
  1220. Default is False.
  1221. check_finite : bool, optional
  1222. Whether to check that the input matrix contains only finite numbers.
  1223. Disabling may give a performance gain, but may result in problems
  1224. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  1225. Returns
  1226. -------
  1227. H : (M, M) ndarray
  1228. Hessenberg form of `a`.
  1229. Q : (M, M) ndarray
  1230. Unitary/orthogonal similarity transformation matrix ``A = Q H Q^H``.
  1231. Only returned if ``calc_q=True``.
  1232. Examples
  1233. --------
  1234. >>> import numpy as np
  1235. >>> from scipy.linalg import hessenberg
  1236. >>> A = np.array([[2, 5, 8, 7], [5, 2, 2, 8], [7, 5, 6, 6], [5, 4, 4, 8]])
  1237. >>> H, Q = hessenberg(A, calc_q=True)
  1238. >>> H
  1239. array([[ 2. , -11.65843866, 1.42005301, 0.25349066],
  1240. [ -9.94987437, 14.53535354, -5.31022304, 2.43081618],
  1241. [ 0. , -1.83299243, 0.38969961, -0.51527034],
  1242. [ 0. , 0. , -3.83189513, 1.07494686]])
  1243. >>> np.allclose(Q @ H @ Q.conj().T - A, np.zeros((4, 4)))
  1244. True
  1245. """
  1246. a1 = _asarray_validated(a, check_finite=check_finite)
  1247. if len(a1.shape) != 2 or (a1.shape[0] != a1.shape[1]):
  1248. raise ValueError('expected square matrix')
  1249. overwrite_a = overwrite_a or (_datacopied(a1, a))
  1250. if a1.size == 0:
  1251. h3 = hessenberg(np.eye(3, dtype=a1.dtype))
  1252. h = np.empty(a1.shape, dtype=h3.dtype)
  1253. if not calc_q:
  1254. return h
  1255. else:
  1256. h3, q3 = hessenberg(np.eye(3, dtype=a1.dtype), calc_q=True)
  1257. q = np.empty(a1.shape, dtype=q3.dtype)
  1258. h = np.empty(a1.shape, dtype=h3.dtype)
  1259. return h, q
  1260. # if 2x2 or smaller: already in Hessenberg
  1261. if a1.shape[0] <= 2:
  1262. if calc_q:
  1263. return a1, eye(a1.shape[0])
  1264. return a1
  1265. gehrd, gebal, gehrd_lwork = get_lapack_funcs(('gehrd', 'gebal',
  1266. 'gehrd_lwork'), (a1,))
  1267. ba, lo, hi, pivscale, info = gebal(a1, permute=0, overwrite_a=overwrite_a)
  1268. _check_info(info, 'gebal (hessenberg)', positive=False)
  1269. n = len(a1)
  1270. lwork = _compute_lwork(gehrd_lwork, ba.shape[0], lo=lo, hi=hi)
  1271. hq, tau, info = gehrd(ba, lo=lo, hi=hi, lwork=lwork, overwrite_a=1)
  1272. _check_info(info, 'gehrd (hessenberg)', positive=False)
  1273. h = np.triu(hq, -1)
  1274. if not calc_q:
  1275. return h
  1276. # use orghr/unghr to compute q
  1277. orghr, orghr_lwork = get_lapack_funcs(('orghr', 'orghr_lwork'), (a1,))
  1278. lwork = _compute_lwork(orghr_lwork, n, lo=lo, hi=hi)
  1279. q, info = orghr(a=hq, tau=tau, lo=lo, hi=hi, lwork=lwork, overwrite_a=1)
  1280. _check_info(info, 'orghr (hessenberg)', positive=False)
  1281. return h, q
  1282. def cdf2rdf(w, v):
  1283. """
  1284. Complex diagonal form to real diagonal block form.
  1285. Converts complex eigenvalues ``w`` and eigenvectors ``v`` to real
  1286. eigenvalues in a block diagonal form ``wr`` and the associated real
  1287. eigenvectors ``vr``, such that::
  1288. vr @ wr = X @ vr
  1289. continues to hold, where ``X`` is the original array for which ``w`` and
  1290. ``v`` are the eigenvalues and eigenvectors.
  1291. .. versionadded:: 1.1.0
  1292. Array argument(s) of this function may have additional
  1293. "batch" dimensions prepended to the core shape. In this case, the array is treated
  1294. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  1295. Parameters
  1296. ----------
  1297. w : (..., M) array_like
  1298. Complex or real eigenvalues, an array or stack of arrays
  1299. Conjugate pairs must not be interleaved, else the wrong result
  1300. will be produced. So ``[1+1j, 1, 1-1j]`` will give a correct result,
  1301. but ``[1+1j, 2+1j, 1-1j, 2-1j]`` will not.
  1302. v : (..., M, M) array_like
  1303. Complex or real eigenvectors, a square array or stack of square arrays.
  1304. Returns
  1305. -------
  1306. wr : (..., M, M) ndarray
  1307. Real diagonal block form of eigenvalues
  1308. vr : (..., M, M) ndarray
  1309. Real eigenvectors associated with ``wr``
  1310. See Also
  1311. --------
  1312. eig : Eigenvalues and right eigenvectors for non-symmetric arrays
  1313. rsf2csf : Convert real Schur form to complex Schur form
  1314. Notes
  1315. -----
  1316. ``w``, ``v`` must be the eigenstructure for some *real* matrix ``X``.
  1317. For example, obtained by ``w, v = scipy.linalg.eig(X)`` or
  1318. ``w, v = numpy.linalg.eig(X)`` in which case ``X`` can also represent
  1319. stacked arrays.
  1320. .. versionadded:: 1.1.0
  1321. Examples
  1322. --------
  1323. >>> import numpy as np
  1324. >>> X = np.array([[1, 2, 3], [0, 4, 5], [0, -5, 4]])
  1325. >>> X
  1326. array([[ 1, 2, 3],
  1327. [ 0, 4, 5],
  1328. [ 0, -5, 4]])
  1329. >>> from scipy import linalg
  1330. >>> w, v = linalg.eig(X)
  1331. >>> w
  1332. array([ 1.+0.j, 4.+5.j, 4.-5.j])
  1333. >>> v
  1334. array([[ 1.00000+0.j , -0.01906-0.40016j, -0.01906+0.40016j],
  1335. [ 0.00000+0.j , 0.00000-0.64788j, 0.00000+0.64788j],
  1336. [ 0.00000+0.j , 0.64788+0.j , 0.64788-0.j ]])
  1337. >>> wr, vr = linalg.cdf2rdf(w, v)
  1338. >>> wr
  1339. array([[ 1., 0., 0.],
  1340. [ 0., 4., 5.],
  1341. [ 0., -5., 4.]])
  1342. >>> vr
  1343. array([[ 1. , 0.40016, -0.01906],
  1344. [ 0. , 0.64788, 0. ],
  1345. [ 0. , 0. , 0.64788]])
  1346. >>> vr @ wr
  1347. array([[ 1. , 1.69593, 1.9246 ],
  1348. [ 0. , 2.59153, 3.23942],
  1349. [ 0. , -3.23942, 2.59153]])
  1350. >>> X @ vr
  1351. array([[ 1. , 1.69593, 1.9246 ],
  1352. [ 0. , 2.59153, 3.23942],
  1353. [ 0. , -3.23942, 2.59153]])
  1354. """
  1355. w, v = _asarray_validated(w), _asarray_validated(v)
  1356. # check dimensions
  1357. if w.ndim < 1:
  1358. raise ValueError('expected w to be at least 1D')
  1359. if v.ndim < 2:
  1360. raise ValueError('expected v to be at least 2D')
  1361. if v.ndim != w.ndim + 1:
  1362. raise ValueError('expected eigenvectors array to have exactly one '
  1363. 'dimension more than eigenvalues array')
  1364. # check shapes
  1365. n = w.shape[-1]
  1366. M = w.shape[:-1]
  1367. if v.shape[-2] != v.shape[-1]:
  1368. raise ValueError('expected v to be a square matrix or stacked square '
  1369. 'matrices: v.shape[-2] = v.shape[-1]')
  1370. if v.shape[-1] != n:
  1371. raise ValueError('expected the same number of eigenvalues as '
  1372. 'eigenvectors')
  1373. # get indices for each first pair of complex eigenvalues
  1374. complex_mask = iscomplex(w)
  1375. n_complex = complex_mask.sum(axis=-1)
  1376. # check if all complex eigenvalues have conjugate pairs
  1377. if not (n_complex % 2 == 0).all():
  1378. raise ValueError('expected complex-conjugate pairs of eigenvalues')
  1379. # find complex indices
  1380. idx = nonzero(complex_mask)
  1381. idx_stack = idx[:-1]
  1382. idx_elem = idx[-1]
  1383. # filter them to conjugate indices, assuming pairs are not interleaved
  1384. j = idx_elem[0::2]
  1385. k = idx_elem[1::2]
  1386. stack_ind = ()
  1387. for i in idx_stack:
  1388. # should never happen, assuming nonzero orders by the last axis
  1389. assert (i[0::2] == i[1::2]).all(), \
  1390. "Conjugate pair spanned different arrays!"
  1391. stack_ind += (i[0::2],)
  1392. # all eigenvalues to diagonal form
  1393. wr = zeros(M + (n, n), dtype=w.real.dtype)
  1394. di = range(n)
  1395. wr[..., di, di] = w.real
  1396. # complex eigenvalues to real block diagonal form
  1397. wr[stack_ind + (j, k)] = w[stack_ind + (j,)].imag
  1398. wr[stack_ind + (k, j)] = w[stack_ind + (k,)].imag
  1399. # compute real eigenvectors associated with real block diagonal eigenvalues
  1400. u = zeros(M + (n, n), dtype=np.cdouble)
  1401. u[..., di, di] = 1.0
  1402. u[stack_ind + (j, j)] = 0.5j
  1403. u[stack_ind + (j, k)] = 0.5
  1404. u[stack_ind + (k, j)] = -0.5j
  1405. u[stack_ind + (k, k)] = 0.5
  1406. # multiply matrices v and u (equivalent to v @ u)
  1407. vr = einsum('...ij,...jk->...ik', v, u).real
  1408. return wr, vr