decomp.py 52 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431
  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 were
  12. # moved to their own files. Still in this file are functions for eigenstuff
  13. # and for the Hessenberg form.
  14. from __future__ import division, print_function, absolute_import
  15. __all__ = ['eig', 'eigvals', 'eigh', 'eigvalsh',
  16. 'eig_banded', 'eigvals_banded',
  17. 'eigh_tridiagonal', 'eigvalsh_tridiagonal', 'hessenberg', 'cdf2rdf']
  18. import numpy
  19. from numpy import (array, isfinite, inexact, nonzero, iscomplexobj, cast,
  20. flatnonzero, conj, asarray, argsort, empty, newaxis,
  21. argwhere, iscomplex, eye, zeros, einsum)
  22. # Local imports
  23. from scipy._lib.six import xrange
  24. from scipy._lib._util import _asarray_validated
  25. from scipy._lib.six import string_types
  26. from .misc import LinAlgError, _datacopied, norm
  27. from .lapack import get_lapack_funcs, _compute_lwork
  28. _I = cast['F'](1j)
  29. def _make_complex_eigvecs(w, vin, dtype):
  30. """
  31. Produce complex-valued eigenvectors from LAPACK DGGEV real-valued output
  32. """
  33. # - see LAPACK man page DGGEV at ALPHAI
  34. v = numpy.array(vin, dtype=dtype)
  35. m = (w.imag > 0)
  36. m[:-1] |= (w.imag[1:] < 0) # workaround for LAPACK bug, cf. ticket #709
  37. for i in flatnonzero(m):
  38. v.imag[:, i] = vin[:, i+1]
  39. conj(v[:, i], v[:, i+1])
  40. return v
  41. def _make_eigvals(alpha, beta, homogeneous_eigvals):
  42. if homogeneous_eigvals:
  43. if beta is None:
  44. return numpy.vstack((alpha, numpy.ones_like(alpha)))
  45. else:
  46. return numpy.vstack((alpha, beta))
  47. else:
  48. if beta is None:
  49. return alpha
  50. else:
  51. w = numpy.empty_like(alpha)
  52. alpha_zero = (alpha == 0)
  53. beta_zero = (beta == 0)
  54. beta_nonzero = ~beta_zero
  55. w[beta_nonzero] = alpha[beta_nonzero]/beta[beta_nonzero]
  56. # Use numpy.inf for complex values too since
  57. # 1/numpy.inf = 0, i.e. it correctly behaves as projective
  58. # infinity.
  59. w[~alpha_zero & beta_zero] = numpy.inf
  60. if numpy.all(alpha.imag == 0):
  61. w[alpha_zero & beta_zero] = numpy.nan
  62. else:
  63. w[alpha_zero & beta_zero] = complex(numpy.nan, numpy.nan)
  64. return w
  65. def _geneig(a1, b1, left, right, overwrite_a, overwrite_b,
  66. homogeneous_eigvals):
  67. ggev, = get_lapack_funcs(('ggev',), (a1, b1))
  68. cvl, cvr = left, right
  69. res = ggev(a1, b1, lwork=-1)
  70. lwork = res[-2][0].real.astype(numpy.int)
  71. if ggev.typecode in 'cz':
  72. alpha, beta, vl, vr, work, info = ggev(a1, b1, cvl, cvr, lwork,
  73. overwrite_a, overwrite_b)
  74. w = _make_eigvals(alpha, beta, homogeneous_eigvals)
  75. else:
  76. alphar, alphai, beta, vl, vr, work, info = ggev(a1, b1, cvl, cvr,
  77. lwork, overwrite_a,
  78. overwrite_b)
  79. alpha = alphar + _I * alphai
  80. w = _make_eigvals(alpha, beta, homogeneous_eigvals)
  81. _check_info(info, 'generalized eig algorithm (ggev)')
  82. only_real = numpy.all(w.imag == 0.0)
  83. if not (ggev.typecode in 'cz' or only_real):
  84. t = w.dtype.char
  85. if left:
  86. vl = _make_complex_eigvecs(w, vl, t)
  87. if right:
  88. vr = _make_complex_eigvecs(w, vr, t)
  89. # the eigenvectors returned by the lapack function are NOT normalized
  90. for i in xrange(vr.shape[0]):
  91. if right:
  92. vr[:, i] /= norm(vr[:, i])
  93. if left:
  94. vl[:, i] /= norm(vl[:, i])
  95. if not (left or right):
  96. return w
  97. if left:
  98. if right:
  99. return w, vl, vr
  100. return w, vl
  101. return w, vr
  102. def eig(a, b=None, left=False, right=True, overwrite_a=False,
  103. overwrite_b=False, check_finite=True, homogeneous_eigvals=False):
  104. """
  105. Solve an ordinary or generalized eigenvalue problem of a square matrix.
  106. Find eigenvalues w and right or left eigenvectors of a general matrix::
  107. a vr[:,i] = w[i] b vr[:,i]
  108. a.H vl[:,i] = w[i].conj() b.H vl[:,i]
  109. where ``.H`` is the Hermitian conjugation.
  110. Parameters
  111. ----------
  112. a : (M, M) array_like
  113. A complex or real matrix whose eigenvalues and eigenvectors
  114. will be computed.
  115. b : (M, M) array_like, optional
  116. Right-hand side matrix in a generalized eigenvalue problem.
  117. Default is None, identity matrix is assumed.
  118. left : bool, optional
  119. Whether to calculate and return left eigenvectors. Default is False.
  120. right : bool, optional
  121. Whether to calculate and return right eigenvectors. Default is True.
  122. overwrite_a : bool, optional
  123. Whether to overwrite `a`; may improve performance. Default is False.
  124. overwrite_b : bool, optional
  125. Whether to overwrite `b`; may improve performance. Default is False.
  126. check_finite : bool, optional
  127. Whether to check that the input matrices contain only finite numbers.
  128. Disabling may give a performance gain, but may result in problems
  129. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  130. homogeneous_eigvals : bool, optional
  131. If True, return the eigenvalues in homogeneous coordinates.
  132. In this case ``w`` is a (2, M) array so that::
  133. w[1,i] a vr[:,i] = w[0,i] b vr[:,i]
  134. Default is False.
  135. Returns
  136. -------
  137. w : (M,) or (2, M) double or complex ndarray
  138. The eigenvalues, each repeated according to its
  139. multiplicity. The shape is (M,) unless
  140. ``homogeneous_eigvals=True``.
  141. vl : (M, M) double or complex ndarray
  142. The normalized left eigenvector corresponding to the eigenvalue
  143. ``w[i]`` is the column vl[:,i]. Only returned if ``left=True``.
  144. vr : (M, M) double or complex ndarray
  145. The normalized right eigenvector corresponding to the eigenvalue
  146. ``w[i]`` is the column ``vr[:,i]``. Only returned if ``right=True``.
  147. Raises
  148. ------
  149. LinAlgError
  150. If eigenvalue computation does not converge.
  151. See Also
  152. --------
  153. eigvals : eigenvalues of general arrays
  154. eigh : Eigenvalues and right eigenvectors for symmetric/Hermitian arrays.
  155. eig_banded : eigenvalues and right eigenvectors for symmetric/Hermitian
  156. band matrices
  157. eigh_tridiagonal : eigenvalues and right eiegenvectors for
  158. symmetric/Hermitian tridiagonal matrices
  159. Examples
  160. --------
  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. overwrite_a = overwrite_a or (_datacopied(a1, a))
  186. if b is not None:
  187. b1 = _asarray_validated(b, check_finite=check_finite)
  188. overwrite_b = overwrite_b or _datacopied(b1, b)
  189. if len(b1.shape) != 2 or b1.shape[0] != b1.shape[1]:
  190. raise ValueError('expected square matrix')
  191. if b1.shape != a1.shape:
  192. raise ValueError('a and b must have the same shape')
  193. return _geneig(a1, b1, left, right, overwrite_a, overwrite_b,
  194. homogeneous_eigvals)
  195. geev, geev_lwork = get_lapack_funcs(('geev', 'geev_lwork'), (a1,))
  196. compute_vl, compute_vr = left, right
  197. lwork = _compute_lwork(geev_lwork, a1.shape[0],
  198. compute_vl=compute_vl,
  199. compute_vr=compute_vr)
  200. if geev.typecode in 'cz':
  201. w, vl, vr, info = geev(a1, lwork=lwork,
  202. compute_vl=compute_vl,
  203. compute_vr=compute_vr,
  204. overwrite_a=overwrite_a)
  205. w = _make_eigvals(w, None, homogeneous_eigvals)
  206. else:
  207. wr, wi, vl, vr, info = geev(a1, lwork=lwork,
  208. compute_vl=compute_vl,
  209. compute_vr=compute_vr,
  210. overwrite_a=overwrite_a)
  211. t = {'f': 'F', 'd': 'D'}[wr.dtype.char]
  212. w = wr + _I * wi
  213. w = _make_eigvals(w, None, homogeneous_eigvals)
  214. _check_info(info, 'eig algorithm (geev)',
  215. positive='did not converge (only eigenvalues '
  216. 'with order >= %d have converged)')
  217. only_real = numpy.all(w.imag == 0.0)
  218. if not (geev.typecode in 'cz' or only_real):
  219. t = w.dtype.char
  220. if left:
  221. vl = _make_complex_eigvecs(w, vl, t)
  222. if right:
  223. vr = _make_complex_eigvecs(w, vr, t)
  224. if not (left or right):
  225. return w
  226. if left:
  227. if right:
  228. return w, vl, vr
  229. return w, vl
  230. return w, vr
  231. def eigh(a, b=None, lower=True, eigvals_only=False, overwrite_a=False,
  232. overwrite_b=False, turbo=True, eigvals=None, type=1,
  233. check_finite=True):
  234. """
  235. Solve an ordinary or generalized eigenvalue problem for a complex
  236. Hermitian or real symmetric matrix.
  237. Find eigenvalues w and optionally eigenvectors v of matrix `a`, where
  238. `b` is positive definite::
  239. a v[:,i] = w[i] b v[:,i]
  240. v[i,:].conj() a v[:,i] = w[i]
  241. v[i,:].conj() b v[:,i] = 1
  242. Parameters
  243. ----------
  244. a : (M, M) array_like
  245. A complex Hermitian or real symmetric matrix whose eigenvalues and
  246. eigenvectors will be computed.
  247. b : (M, M) array_like, optional
  248. A complex Hermitian or real symmetric definite positive matrix in.
  249. If omitted, identity matrix is assumed.
  250. lower : bool, optional
  251. Whether the pertinent array data is taken from the lower or upper
  252. triangle of `a`. (Default: lower)
  253. eigvals_only : bool, optional
  254. Whether to calculate only eigenvalues and no eigenvectors.
  255. (Default: both are calculated)
  256. turbo : bool, optional
  257. Use divide and conquer algorithm (faster but expensive in memory,
  258. only for generalized eigenvalue problem and if eigvals=None)
  259. eigvals : tuple (lo, hi), optional
  260. Indexes of the smallest and largest (in ascending order) eigenvalues
  261. and corresponding eigenvectors to be returned: 0 <= lo <= hi <= M-1.
  262. If omitted, all eigenvalues and eigenvectors are returned.
  263. type : int, optional
  264. Specifies the problem type to be solved:
  265. type = 1: a v[:,i] = w[i] b v[:,i]
  266. type = 2: a b v[:,i] = w[i] v[:,i]
  267. type = 3: b a v[:,i] = w[i] v[:,i]
  268. overwrite_a : bool, optional
  269. Whether to overwrite data in `a` (may improve performance)
  270. overwrite_b : bool, optional
  271. Whether to overwrite data in `b` (may improve performance)
  272. check_finite : bool, optional
  273. Whether to check that the input matrices contain only finite numbers.
  274. Disabling may give a performance gain, but may result in problems
  275. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  276. Returns
  277. -------
  278. w : (N,) float ndarray
  279. The N (1<=N<=M) selected eigenvalues, in ascending order, each
  280. repeated according to its multiplicity.
  281. v : (M, N) complex ndarray
  282. (if eigvals_only == False)
  283. The normalized selected eigenvector corresponding to the
  284. eigenvalue w[i] is the column v[:,i].
  285. Normalization:
  286. type 1 and 3: v.conj() a v = w
  287. type 2: inv(v).conj() a inv(v) = w
  288. type = 1 or 2: v.conj() b v = I
  289. type = 3: v.conj() inv(b) v = I
  290. Raises
  291. ------
  292. LinAlgError
  293. If eigenvalue computation does not converge,
  294. an error occurred, or b matrix is not definite positive. Note that
  295. if input matrices are not symmetric or hermitian, no error is reported
  296. but results will be wrong.
  297. See Also
  298. --------
  299. eigvalsh : eigenvalues of symmetric or Hermitian arrays
  300. eig : eigenvalues and right eigenvectors for non-symmetric arrays
  301. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  302. eigh_tridiagonal : eigenvalues and right eiegenvectors for
  303. symmetric/Hermitian tridiagonal matrices
  304. Notes
  305. -----
  306. This function does not check the input array for being hermitian/symmetric
  307. in order to allow for representing arrays with only their upper/lower
  308. triangular parts.
  309. Examples
  310. --------
  311. >>> from scipy.linalg import eigh
  312. >>> A = np.array([[6, 3, 1, 5], [3, 0, 5, 1], [1, 5, 6, 2], [5, 1, 2, 2]])
  313. >>> w, v = eigh(A)
  314. >>> np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4)))
  315. True
  316. """
  317. a1 = _asarray_validated(a, check_finite=check_finite)
  318. if len(a1.shape) != 2 or a1.shape[0] != a1.shape[1]:
  319. raise ValueError('expected square matrix')
  320. overwrite_a = overwrite_a or (_datacopied(a1, a))
  321. if iscomplexobj(a1):
  322. cplx = True
  323. else:
  324. cplx = False
  325. if b is not None:
  326. b1 = _asarray_validated(b, check_finite=check_finite)
  327. overwrite_b = overwrite_b or _datacopied(b1, b)
  328. if len(b1.shape) != 2 or b1.shape[0] != b1.shape[1]:
  329. raise ValueError('expected square matrix')
  330. if b1.shape != a1.shape:
  331. raise ValueError("wrong b dimensions %s, should "
  332. "be %s" % (str(b1.shape), str(a1.shape)))
  333. if iscomplexobj(b1):
  334. cplx = True
  335. else:
  336. cplx = cplx or False
  337. else:
  338. b1 = None
  339. # Set job for fortran routines
  340. _job = (eigvals_only and 'N') or 'V'
  341. # port eigenvalue range from python to fortran convention
  342. if eigvals is not None:
  343. lo, hi = eigvals
  344. if lo < 0 or hi >= a1.shape[0]:
  345. raise ValueError('The eigenvalue range specified is not valid.\n'
  346. 'Valid range is [%s,%s]' % (0, a1.shape[0]-1))
  347. lo += 1
  348. hi += 1
  349. eigvals = (lo, hi)
  350. # set lower
  351. if lower:
  352. uplo = 'L'
  353. else:
  354. uplo = 'U'
  355. # fix prefix for lapack routines
  356. if cplx:
  357. pfx = 'he'
  358. else:
  359. pfx = 'sy'
  360. # Standard Eigenvalue Problem
  361. # Use '*evr' routines
  362. # FIXME: implement calculation of optimal lwork
  363. # for all lapack routines
  364. if b1 is None:
  365. driver = pfx+'evr'
  366. (evr,) = get_lapack_funcs((driver,), (a1,))
  367. if eigvals is None:
  368. w, v, info = evr(a1, uplo=uplo, jobz=_job, range="A", il=1,
  369. iu=a1.shape[0], overwrite_a=overwrite_a)
  370. else:
  371. (lo, hi) = eigvals
  372. w_tot, v, info = evr(a1, uplo=uplo, jobz=_job, range="I",
  373. il=lo, iu=hi, overwrite_a=overwrite_a)
  374. w = w_tot[0:hi-lo+1]
  375. # Generalized Eigenvalue Problem
  376. else:
  377. # Use '*gvx' routines if range is specified
  378. if eigvals is not None:
  379. driver = pfx+'gvx'
  380. (gvx,) = get_lapack_funcs((driver,), (a1, b1))
  381. (lo, hi) = eigvals
  382. w_tot, v, ifail, info = gvx(a1, b1, uplo=uplo, iu=hi,
  383. itype=type, jobz=_job, il=lo,
  384. overwrite_a=overwrite_a,
  385. overwrite_b=overwrite_b)
  386. w = w_tot[0:hi-lo+1]
  387. # Use '*gvd' routine if turbo is on and no eigvals are specified
  388. elif turbo:
  389. driver = pfx+'gvd'
  390. (gvd,) = get_lapack_funcs((driver,), (a1, b1))
  391. v, w, info = gvd(a1, b1, uplo=uplo, itype=type, jobz=_job,
  392. overwrite_a=overwrite_a,
  393. overwrite_b=overwrite_b)
  394. # Use '*gv' routine if turbo is off and no eigvals are specified
  395. else:
  396. driver = pfx+'gv'
  397. (gv,) = get_lapack_funcs((driver,), (a1, b1))
  398. v, w, info = gv(a1, b1, uplo=uplo, itype=type, jobz=_job,
  399. overwrite_a=overwrite_a,
  400. overwrite_b=overwrite_b)
  401. # Check if we had a successful exit
  402. if info == 0:
  403. if eigvals_only:
  404. return w
  405. else:
  406. return w, v
  407. _check_info(info, driver, positive=False) # triage more specifically
  408. if info > 0 and b1 is None:
  409. raise LinAlgError("unrecoverable internal error.")
  410. # The algorithm failed to converge.
  411. elif 0 < info <= b1.shape[0]:
  412. if eigvals is not None:
  413. raise LinAlgError("the eigenvectors %s failed to"
  414. " converge." % nonzero(ifail)-1)
  415. else:
  416. raise LinAlgError("internal fortran routine failed to converge: "
  417. "%i off-diagonal elements of an "
  418. "intermediate tridiagonal form did not converge"
  419. " to zero." % info)
  420. # This occurs when b is not positive definite
  421. else:
  422. raise LinAlgError("the leading minor of order %i"
  423. " of 'b' is not positive definite. The"
  424. " factorization of 'b' could not be completed"
  425. " and no eigenvalues or eigenvectors were"
  426. " computed." % (info-b1.shape[0]))
  427. _conv_dict = {0: 0, 1: 1, 2: 2,
  428. 'all': 0, 'value': 1, 'index': 2,
  429. 'a': 0, 'v': 1, 'i': 2}
  430. def _check_select(select, select_range, max_ev, max_len):
  431. """Check that select is valid, convert to Fortran style."""
  432. if isinstance(select, string_types):
  433. select = select.lower()
  434. try:
  435. select = _conv_dict[select]
  436. except KeyError:
  437. raise ValueError('invalid argument for select')
  438. vl, vu = 0., 1.
  439. il = iu = 1
  440. if select != 0: # (non-all)
  441. sr = asarray(select_range)
  442. if sr.ndim != 1 or sr.size != 2 or sr[1] < sr[0]:
  443. raise ValueError('select_range must be a 2-element array-like '
  444. 'in nondecreasing order')
  445. if select == 1: # (value)
  446. vl, vu = sr
  447. if max_ev == 0:
  448. max_ev = max_len
  449. else: # 2 (index)
  450. if sr.dtype.char.lower() not in 'hilqp':
  451. raise ValueError('when using select="i", select_range must '
  452. 'contain integers, got dtype %s (%s)'
  453. % (sr.dtype, sr.dtype.char))
  454. # translate Python (0 ... N-1) into Fortran (1 ... N) with + 1
  455. il, iu = sr + 1
  456. if min(il, iu) < 1 or max(il, iu) > max_len:
  457. raise ValueError('select_range out of bounds')
  458. max_ev = iu - il + 1
  459. return select, vl, vu, il, iu, max_ev
  460. def eig_banded(a_band, lower=False, eigvals_only=False, overwrite_a_band=False,
  461. select='a', select_range=None, max_ev=0, check_finite=True):
  462. """
  463. Solve real symmetric or complex hermitian band matrix eigenvalue problem.
  464. Find eigenvalues w and optionally right eigenvectors v of a::
  465. a v[:,i] = w[i] v[:,i]
  466. v.H v = identity
  467. The matrix a is stored in a_band either in lower diagonal or upper
  468. diagonal ordered form:
  469. a_band[u + i - j, j] == a[i,j] (if upper form; i <= j)
  470. a_band[ i - j, j] == a[i,j] (if lower form; i >= j)
  471. where u is the number of bands above the diagonal.
  472. Example of a_band (shape of a is (6,6), u=2)::
  473. upper form:
  474. * * a02 a13 a24 a35
  475. * a01 a12 a23 a34 a45
  476. a00 a11 a22 a33 a44 a55
  477. lower form:
  478. a00 a11 a22 a33 a44 a55
  479. a10 a21 a32 a43 a54 *
  480. a20 a31 a42 a53 * *
  481. Cells marked with * are not used.
  482. Parameters
  483. ----------
  484. a_band : (u+1, M) array_like
  485. The bands of the M by M matrix a.
  486. lower : bool, optional
  487. Is the matrix in the lower form. (Default is upper form)
  488. eigvals_only : bool, optional
  489. Compute only the eigenvalues and no eigenvectors.
  490. (Default: calculate also eigenvectors)
  491. overwrite_a_band : bool, optional
  492. Discard data in a_band (may enhance performance)
  493. select : {'a', 'v', 'i'}, optional
  494. Which eigenvalues to calculate
  495. ====== ========================================
  496. select calculated
  497. ====== ========================================
  498. 'a' All eigenvalues
  499. 'v' Eigenvalues in the interval (min, max]
  500. 'i' Eigenvalues with indices min <= i <= max
  501. ====== ========================================
  502. select_range : (min, max), optional
  503. Range of selected eigenvalues
  504. max_ev : int, optional
  505. For select=='v', maximum number of eigenvalues expected.
  506. For other values of select, has no meaning.
  507. In doubt, leave this parameter untouched.
  508. check_finite : bool, optional
  509. Whether to check that the input matrix contains only finite numbers.
  510. Disabling may give a performance gain, but may result in problems
  511. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  512. Returns
  513. -------
  514. w : (M,) ndarray
  515. The eigenvalues, in ascending order, each repeated according to its
  516. multiplicity.
  517. v : (M, M) float or complex ndarray
  518. The normalized eigenvector corresponding to the eigenvalue w[i] is
  519. the column v[:,i].
  520. Raises
  521. ------
  522. LinAlgError
  523. If eigenvalue computation does not converge.
  524. See Also
  525. --------
  526. eigvals_banded : eigenvalues for symmetric/Hermitian band matrices
  527. eig : eigenvalues and right eigenvectors of general arrays.
  528. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  529. eigh_tridiagonal : eigenvalues and right eiegenvectors for
  530. symmetric/Hermitian tridiagonal matrices
  531. Examples
  532. --------
  533. >>> from scipy.linalg import eig_banded
  534. >>> A = np.array([[1, 5, 2, 0], [5, 2, 5, 2], [2, 5, 3, 5], [0, 2, 5, 4]])
  535. >>> Ab = np.array([[1, 2, 3, 4], [5, 5, 5, 0], [2, 2, 0, 0]])
  536. >>> w, v = eig_banded(Ab, lower=True)
  537. >>> np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4)))
  538. True
  539. >>> w = eig_banded(Ab, lower=True, eigvals_only=True)
  540. >>> w
  541. array([-4.26200532, -2.22987175, 3.95222349, 12.53965359])
  542. Request only the eigenvalues between ``[-3, 4]``
  543. >>> w, v = eig_banded(Ab, lower=True, select='v', select_range=[-3, 4])
  544. >>> w
  545. array([-2.22987175, 3.95222349])
  546. """
  547. if eigvals_only or overwrite_a_band:
  548. a1 = _asarray_validated(a_band, check_finite=check_finite)
  549. overwrite_a_band = overwrite_a_band or (_datacopied(a1, a_band))
  550. else:
  551. a1 = array(a_band)
  552. if issubclass(a1.dtype.type, inexact) and not isfinite(a1).all():
  553. raise ValueError("array must not contain infs or NaNs")
  554. overwrite_a_band = 1
  555. if len(a1.shape) != 2:
  556. raise ValueError('expected two-dimensional array')
  557. select, vl, vu, il, iu, max_ev = _check_select(
  558. select, select_range, max_ev, a1.shape[1])
  559. del select_range
  560. if select == 0:
  561. if a1.dtype.char in 'GFD':
  562. # FIXME: implement this somewhen, for now go with builtin values
  563. # FIXME: calc optimal lwork by calling ?hbevd(lwork=-1)
  564. # or by using calc_lwork.f ???
  565. # lwork = calc_lwork.hbevd(bevd.typecode, a1.shape[0], lower)
  566. internal_name = 'hbevd'
  567. else: # a1.dtype.char in 'fd':
  568. # FIXME: implement this somewhen, for now go with builtin values
  569. # see above
  570. # lwork = calc_lwork.sbevd(bevd.typecode, a1.shape[0], lower)
  571. internal_name = 'sbevd'
  572. bevd, = get_lapack_funcs((internal_name,), (a1,))
  573. w, v, info = bevd(a1, compute_v=not eigvals_only,
  574. lower=lower, overwrite_ab=overwrite_a_band)
  575. else: # select in [1, 2]
  576. if eigvals_only:
  577. max_ev = 1
  578. # calculate optimal abstol for dsbevx (see manpage)
  579. if a1.dtype.char in 'fF': # single precision
  580. lamch, = get_lapack_funcs(('lamch',), (array(0, dtype='f'),))
  581. else:
  582. lamch, = get_lapack_funcs(('lamch',), (array(0, dtype='d'),))
  583. abstol = 2 * lamch('s')
  584. if a1.dtype.char in 'GFD':
  585. internal_name = 'hbevx'
  586. else: # a1.dtype.char in 'gfd'
  587. internal_name = 'sbevx'
  588. bevx, = get_lapack_funcs((internal_name,), (a1,))
  589. w, v, m, ifail, info = bevx(
  590. a1, vl, vu, il, iu, compute_v=not eigvals_only, mmax=max_ev,
  591. range=select, lower=lower, overwrite_ab=overwrite_a_band,
  592. abstol=abstol)
  593. # crop off w and v
  594. w = w[:m]
  595. if not eigvals_only:
  596. v = v[:, :m]
  597. _check_info(info, internal_name)
  598. if eigvals_only:
  599. return w
  600. return w, v
  601. def eigvals(a, b=None, overwrite_a=False, check_finite=True,
  602. homogeneous_eigvals=False):
  603. """
  604. Compute eigenvalues from an ordinary or generalized eigenvalue problem.
  605. Find eigenvalues of a general matrix::
  606. a vr[:,i] = w[i] b vr[:,i]
  607. Parameters
  608. ----------
  609. a : (M, M) array_like
  610. A complex or real matrix whose eigenvalues and eigenvectors
  611. will be computed.
  612. b : (M, M) array_like, optional
  613. Right-hand side matrix in a generalized eigenvalue problem.
  614. If omitted, identity matrix is assumed.
  615. overwrite_a : bool, optional
  616. Whether to overwrite data in a (may improve performance)
  617. check_finite : bool, optional
  618. Whether to check that the input matrices contain only finite numbers.
  619. Disabling may give a performance gain, but may result in problems
  620. (crashes, non-termination) if the inputs do contain infinities
  621. or NaNs.
  622. homogeneous_eigvals : bool, optional
  623. If True, return the eigenvalues in homogeneous coordinates.
  624. In this case ``w`` is a (2, M) array so that::
  625. w[1,i] a vr[:,i] = w[0,i] b vr[:,i]
  626. Default is False.
  627. Returns
  628. -------
  629. w : (M,) or (2, M) double or complex ndarray
  630. The eigenvalues, each repeated according to its multiplicity
  631. but not in any specific order. The shape is (M,) unless
  632. ``homogeneous_eigvals=True``.
  633. Raises
  634. ------
  635. LinAlgError
  636. If eigenvalue computation does not converge
  637. See Also
  638. --------
  639. eig : eigenvalues and right eigenvectors of general arrays.
  640. eigvalsh : eigenvalues of symmetric or Hermitian arrays
  641. eigvals_banded : eigenvalues for symmetric/Hermitian band matrices
  642. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  643. matrices
  644. Examples
  645. --------
  646. >>> from scipy import linalg
  647. >>> a = np.array([[0., -1.], [1., 0.]])
  648. >>> linalg.eigvals(a)
  649. array([0.+1.j, 0.-1.j])
  650. >>> b = np.array([[0., 1.], [1., 1.]])
  651. >>> linalg.eigvals(a, b)
  652. array([ 1.+0.j, -1.+0.j])
  653. >>> a = np.array([[3., 0., 0.], [0., 8., 0.], [0., 0., 7.]])
  654. >>> linalg.eigvals(a, homogeneous_eigvals=True)
  655. array([[3.+0.j, 8.+0.j, 7.+0.j],
  656. [1.+0.j, 1.+0.j, 1.+0.j]])
  657. """
  658. return eig(a, b=b, left=0, right=0, overwrite_a=overwrite_a,
  659. check_finite=check_finite,
  660. homogeneous_eigvals=homogeneous_eigvals)
  661. def eigvalsh(a, b=None, lower=True, overwrite_a=False,
  662. overwrite_b=False, turbo=True, eigvals=None, type=1,
  663. check_finite=True):
  664. """
  665. Solve an ordinary or generalized eigenvalue problem for a complex
  666. Hermitian or real symmetric matrix.
  667. Find eigenvalues w of matrix a, where b is positive definite::
  668. a v[:,i] = w[i] b v[:,i]
  669. v[i,:].conj() a v[:,i] = w[i]
  670. v[i,:].conj() b v[:,i] = 1
  671. Parameters
  672. ----------
  673. a : (M, M) array_like
  674. A complex Hermitian or real symmetric matrix whose eigenvalues and
  675. eigenvectors will be computed.
  676. b : (M, M) array_like, optional
  677. A complex Hermitian or real symmetric definite positive matrix in.
  678. If omitted, identity matrix is assumed.
  679. lower : bool, optional
  680. Whether the pertinent array data is taken from the lower or upper
  681. triangle of `a`. (Default: lower)
  682. turbo : bool, optional
  683. Use divide and conquer algorithm (faster but expensive in memory,
  684. only for generalized eigenvalue problem and if eigvals=None)
  685. eigvals : tuple (lo, hi), optional
  686. Indexes of the smallest and largest (in ascending order) eigenvalues
  687. and corresponding eigenvectors to be returned: 0 <= lo < hi <= M-1.
  688. If omitted, all eigenvalues and eigenvectors are returned.
  689. type : int, optional
  690. Specifies the problem type to be solved:
  691. type = 1: a v[:,i] = w[i] b v[:,i]
  692. type = 2: a b v[:,i] = w[i] v[:,i]
  693. type = 3: b a v[:,i] = w[i] v[:,i]
  694. overwrite_a : bool, optional
  695. Whether to overwrite data in `a` (may improve performance)
  696. overwrite_b : bool, optional
  697. Whether to overwrite data in `b` (may improve performance)
  698. check_finite : bool, optional
  699. Whether to check that the input matrices contain only finite numbers.
  700. Disabling may give a performance gain, but may result in problems
  701. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  702. Returns
  703. -------
  704. w : (N,) float ndarray
  705. The N (1<=N<=M) selected eigenvalues, in ascending order, each
  706. repeated according to its multiplicity.
  707. Raises
  708. ------
  709. LinAlgError
  710. If eigenvalue computation does not converge,
  711. an error occurred, or b matrix is not definite positive. Note that
  712. if input matrices are not symmetric or hermitian, no error is reported
  713. but results will be wrong.
  714. See Also
  715. --------
  716. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  717. eigvals : eigenvalues of general arrays
  718. eigvals_banded : eigenvalues for symmetric/Hermitian band matrices
  719. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  720. matrices
  721. Notes
  722. -----
  723. This function does not check the input array for being hermitian/symmetric
  724. in order to allow for representing arrays with only their upper/lower
  725. triangular parts.
  726. Examples
  727. --------
  728. >>> from scipy.linalg import eigvalsh
  729. >>> A = np.array([[6, 3, 1, 5], [3, 0, 5, 1], [1, 5, 6, 2], [5, 1, 2, 2]])
  730. >>> w = eigvalsh(A)
  731. >>> w
  732. array([-3.74637491, -0.76263923, 6.08502336, 12.42399079])
  733. """
  734. return eigh(a, b=b, lower=lower, eigvals_only=True,
  735. overwrite_a=overwrite_a, overwrite_b=overwrite_b,
  736. turbo=turbo, eigvals=eigvals, type=type,
  737. check_finite=check_finite)
  738. def eigvals_banded(a_band, lower=False, overwrite_a_band=False,
  739. select='a', select_range=None, check_finite=True):
  740. """
  741. Solve real symmetric or complex hermitian band matrix eigenvalue problem.
  742. Find eigenvalues w of a::
  743. a v[:,i] = w[i] v[:,i]
  744. v.H v = identity
  745. The matrix a is stored in a_band either in lower diagonal or upper
  746. diagonal ordered form:
  747. a_band[u + i - j, j] == a[i,j] (if upper form; i <= j)
  748. a_band[ i - j, j] == a[i,j] (if lower form; i >= j)
  749. where u is the number of bands above the diagonal.
  750. Example of a_band (shape of a is (6,6), u=2)::
  751. upper form:
  752. * * a02 a13 a24 a35
  753. * a01 a12 a23 a34 a45
  754. a00 a11 a22 a33 a44 a55
  755. lower form:
  756. a00 a11 a22 a33 a44 a55
  757. a10 a21 a32 a43 a54 *
  758. a20 a31 a42 a53 * *
  759. Cells marked with * are not used.
  760. Parameters
  761. ----------
  762. a_band : (u+1, M) array_like
  763. The bands of the M by M matrix a.
  764. lower : bool, optional
  765. Is the matrix in the lower form. (Default is upper form)
  766. overwrite_a_band : bool, optional
  767. Discard data in a_band (may enhance performance)
  768. select : {'a', 'v', 'i'}, optional
  769. Which eigenvalues to calculate
  770. ====== ========================================
  771. select calculated
  772. ====== ========================================
  773. 'a' All eigenvalues
  774. 'v' Eigenvalues in the interval (min, max]
  775. 'i' Eigenvalues with indices min <= i <= max
  776. ====== ========================================
  777. select_range : (min, max), optional
  778. Range of selected eigenvalues
  779. check_finite : bool, optional
  780. Whether to check that the input matrix contains only finite numbers.
  781. Disabling may give a performance gain, but may result in problems
  782. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  783. Returns
  784. -------
  785. w : (M,) ndarray
  786. The eigenvalues, in ascending order, each repeated according to its
  787. multiplicity.
  788. Raises
  789. ------
  790. LinAlgError
  791. If eigenvalue computation does not converge.
  792. See Also
  793. --------
  794. eig_banded : eigenvalues and right eigenvectors for symmetric/Hermitian
  795. band matrices
  796. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  797. matrices
  798. eigvals : eigenvalues of general arrays
  799. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  800. eig : eigenvalues and right eigenvectors for non-symmetric arrays
  801. Examples
  802. --------
  803. >>> from scipy.linalg import eigvals_banded
  804. >>> A = np.array([[1, 5, 2, 0], [5, 2, 5, 2], [2, 5, 3, 5], [0, 2, 5, 4]])
  805. >>> Ab = np.array([[1, 2, 3, 4], [5, 5, 5, 0], [2, 2, 0, 0]])
  806. >>> w = eigvals_banded(Ab, lower=True)
  807. >>> w
  808. array([-4.26200532, -2.22987175, 3.95222349, 12.53965359])
  809. """
  810. return eig_banded(a_band, lower=lower, eigvals_only=1,
  811. overwrite_a_band=overwrite_a_band, select=select,
  812. select_range=select_range, check_finite=check_finite)
  813. def eigvalsh_tridiagonal(d, e, select='a', select_range=None,
  814. check_finite=True, tol=0., lapack_driver='auto'):
  815. """
  816. Solve eigenvalue problem for a real symmetric tridiagonal matrix.
  817. Find eigenvalues `w` of ``a``::
  818. a v[:,i] = w[i] v[:,i]
  819. v.H v = identity
  820. For a real symmetric matrix ``a`` with diagonal elements `d` and
  821. off-diagonal elements `e`.
  822. Parameters
  823. ----------
  824. d : ndarray, shape (ndim,)
  825. The diagonal elements of the array.
  826. e : ndarray, shape (ndim-1,)
  827. The off-diagonal elements of the array.
  828. select : {'a', 'v', 'i'}, optional
  829. Which eigenvalues to calculate
  830. ====== ========================================
  831. select calculated
  832. ====== ========================================
  833. 'a' All eigenvalues
  834. 'v' Eigenvalues in the interval (min, max]
  835. 'i' Eigenvalues with indices min <= i <= max
  836. ====== ========================================
  837. select_range : (min, max), optional
  838. Range of selected eigenvalues
  839. check_finite : bool, optional
  840. Whether to check that the input matrix contains only finite numbers.
  841. Disabling may give a performance gain, but may result in problems
  842. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  843. tol : float
  844. The absolute tolerance to which each eigenvalue is required
  845. (only used when ``lapack_driver='stebz'``).
  846. An eigenvalue (or cluster) is considered to have converged if it
  847. lies in an interval of this width. If <= 0. (default),
  848. the value ``eps*|a|`` is used where eps is the machine precision,
  849. and ``|a|`` is the 1-norm of the matrix ``a``.
  850. lapack_driver : str
  851. LAPACK function to use, can be 'auto', 'stemr', 'stebz', 'sterf',
  852. or 'stev'. When 'auto' (default), it will use 'stemr' if ``select='a'``
  853. and 'stebz' otherwise. 'sterf' and 'stev' can only be used when
  854. ``select='a'``.
  855. Returns
  856. -------
  857. w : (M,) ndarray
  858. The eigenvalues, in ascending order, each repeated according to its
  859. multiplicity.
  860. Raises
  861. ------
  862. LinAlgError
  863. If eigenvalue computation does not converge.
  864. See Also
  865. --------
  866. eigh_tridiagonal : eigenvalues and right eiegenvectors for
  867. symmetric/Hermitian tridiagonal matrices
  868. Examples
  869. --------
  870. >>> from scipy.linalg import eigvalsh_tridiagonal, eigvalsh
  871. >>> d = 3*np.ones(4)
  872. >>> e = -1*np.ones(3)
  873. >>> w = eigvalsh_tridiagonal(d, e)
  874. >>> A = np.diag(d) + np.diag(e, k=1) + np.diag(e, k=-1)
  875. >>> w2 = eigvalsh(A) # Verify with other eigenvalue routines
  876. >>> np.allclose(w - w2, np.zeros(4))
  877. True
  878. """
  879. return eigh_tridiagonal(
  880. d, e, eigvals_only=True, select=select, select_range=select_range,
  881. check_finite=check_finite, tol=tol, lapack_driver=lapack_driver)
  882. def eigh_tridiagonal(d, e, eigvals_only=False, select='a', select_range=None,
  883. check_finite=True, tol=0., lapack_driver='auto'):
  884. """
  885. Solve eigenvalue problem for a real symmetric tridiagonal matrix.
  886. Find eigenvalues `w` and optionally right eigenvectors `v` of ``a``::
  887. a v[:,i] = w[i] v[:,i]
  888. v.H v = identity
  889. For a real symmetric matrix ``a`` with diagonal elements `d` and
  890. off-diagonal elements `e`.
  891. Parameters
  892. ----------
  893. d : ndarray, shape (ndim,)
  894. The diagonal elements of the array.
  895. e : ndarray, shape (ndim-1,)
  896. The off-diagonal elements of the array.
  897. select : {'a', 'v', 'i'}, optional
  898. Which eigenvalues to calculate
  899. ====== ========================================
  900. select calculated
  901. ====== ========================================
  902. 'a' All eigenvalues
  903. 'v' Eigenvalues in the interval (min, max]
  904. 'i' Eigenvalues with indices min <= i <= max
  905. ====== ========================================
  906. select_range : (min, max), optional
  907. Range of selected eigenvalues
  908. check_finite : bool, optional
  909. Whether to check that the input matrix contains only finite numbers.
  910. Disabling may give a performance gain, but may result in problems
  911. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  912. tol : float
  913. The absolute tolerance to which each eigenvalue is required
  914. (only used when 'stebz' is the `lapack_driver`).
  915. An eigenvalue (or cluster) is considered to have converged if it
  916. lies in an interval of this width. If <= 0. (default),
  917. the value ``eps*|a|`` is used where eps is the machine precision,
  918. and ``|a|`` is the 1-norm of the matrix ``a``.
  919. lapack_driver : str
  920. LAPACK function to use, can be 'auto', 'stemr', 'stebz', 'sterf',
  921. or 'stev'. When 'auto' (default), it will use 'stemr' if ``select='a'``
  922. and 'stebz' otherwise. When 'stebz' is used to find the eigenvalues and
  923. ``eigvals_only=False``, then a second LAPACK call (to ``?STEIN``) is
  924. used to find the corresponding eigenvectors. 'sterf' can only be
  925. used when ``eigvals_only=True`` and ``select='a'``. 'stev' can only
  926. be used when ``select='a'``.
  927. Returns
  928. -------
  929. w : (M,) ndarray
  930. The eigenvalues, in ascending order, each repeated according to its
  931. multiplicity.
  932. v : (M, M) ndarray
  933. The normalized eigenvector corresponding to the eigenvalue ``w[i]`` is
  934. the column ``v[:,i]``.
  935. Raises
  936. ------
  937. LinAlgError
  938. If eigenvalue computation does not converge.
  939. See Also
  940. --------
  941. eigvalsh_tridiagonal : eigenvalues of symmetric/Hermitian tridiagonal
  942. matrices
  943. eig : eigenvalues and right eigenvectors for non-symmetric arrays
  944. eigh : eigenvalues and right eigenvectors for symmetric/Hermitian arrays
  945. eig_banded : eigenvalues and right eigenvectors for symmetric/Hermitian
  946. band matrices
  947. Notes
  948. -----
  949. This function makes use of LAPACK ``S/DSTEMR`` routines.
  950. Examples
  951. --------
  952. >>> from scipy.linalg import eigh_tridiagonal
  953. >>> d = 3*np.ones(4)
  954. >>> e = -1*np.ones(3)
  955. >>> w, v = eigh_tridiagonal(d, e)
  956. >>> A = np.diag(d) + np.diag(e, k=1) + np.diag(e, k=-1)
  957. >>> np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4)))
  958. True
  959. """
  960. d = _asarray_validated(d, check_finite=check_finite)
  961. e = _asarray_validated(e, check_finite=check_finite)
  962. for check in (d, e):
  963. if check.ndim != 1:
  964. raise ValueError('expected one-dimensional array')
  965. if check.dtype.char in 'GFD': # complex
  966. raise TypeError('Only real arrays currently supported')
  967. if d.size != e.size + 1:
  968. raise ValueError('d (%s) must have one more element than e (%s)'
  969. % (d.size, e.size))
  970. select, vl, vu, il, iu, _ = _check_select(
  971. select, select_range, 0, d.size)
  972. if not isinstance(lapack_driver, string_types):
  973. raise TypeError('lapack_driver must be str')
  974. drivers = ('auto', 'stemr', 'sterf', 'stebz', 'stev')
  975. if lapack_driver not in drivers:
  976. raise ValueError('lapack_driver must be one of %s, got %s'
  977. % (drivers, lapack_driver))
  978. if lapack_driver == 'auto':
  979. lapack_driver = 'stemr' if select == 0 else 'stebz'
  980. func, = get_lapack_funcs((lapack_driver,), (d, e))
  981. compute_v = not eigvals_only
  982. if lapack_driver == 'sterf':
  983. if select != 0:
  984. raise ValueError('sterf can only be used when select == "a"')
  985. if not eigvals_only:
  986. raise ValueError('sterf can only be used when eigvals_only is '
  987. 'True')
  988. w, info = func(d, e)
  989. m = len(w)
  990. elif lapack_driver == 'stev':
  991. if select != 0:
  992. raise ValueError('stev can only be used when select == "a"')
  993. w, v, info = func(d, e, compute_v=compute_v)
  994. m = len(w)
  995. elif lapack_driver == 'stebz':
  996. tol = float(tol)
  997. internal_name = 'stebz'
  998. stebz, = get_lapack_funcs((internal_name,), (d, e))
  999. # If getting eigenvectors, needs to be block-ordered (B) instead of
  1000. # matirx-ordered (E), and we will reorder later
  1001. order = 'E' if eigvals_only else 'B'
  1002. m, w, iblock, isplit, info = stebz(d, e, select, vl, vu, il, iu, tol,
  1003. order)
  1004. else: # 'stemr'
  1005. # ?STEMR annoyingly requires size N instead of N-1
  1006. e_ = empty(e.size+1, e.dtype)
  1007. e_[:-1] = e
  1008. stemr_lwork, = get_lapack_funcs(('stemr_lwork',), (d, e))
  1009. lwork, liwork, info = stemr_lwork(d, e_, select, vl, vu, il, iu,
  1010. compute_v=compute_v)
  1011. _check_info(info, 'stemr_lwork')
  1012. m, w, v, info = func(d, e_, select, vl, vu, il, iu,
  1013. compute_v=compute_v, lwork=lwork, liwork=liwork)
  1014. _check_info(info, lapack_driver + ' (eigh_tridiagonal)')
  1015. w = w[:m]
  1016. if eigvals_only:
  1017. return w
  1018. else:
  1019. # Do we still need to compute the eigenvalues?
  1020. if lapack_driver == 'stebz':
  1021. func, = get_lapack_funcs(('stein',), (d, e))
  1022. v, info = func(d, e, w, iblock, isplit)
  1023. _check_info(info, 'stein (eigh_tridiagonal)',
  1024. positive='%d eigenvectors failed to converge')
  1025. # Convert block-order to matrix-order
  1026. order = argsort(w)
  1027. w, v = w[order], v[:, order]
  1028. else:
  1029. v = v[:, :m]
  1030. return w, v
  1031. def _check_info(info, driver, positive='did not converge (LAPACK info=%d)'):
  1032. """Check info return value."""
  1033. if info < 0:
  1034. raise ValueError('illegal value in argument %d of internal %s'
  1035. % (-info, driver))
  1036. if info > 0 and positive:
  1037. raise LinAlgError(("%s " + positive) % (driver, info,))
  1038. def hessenberg(a, calc_q=False, overwrite_a=False, check_finite=True):
  1039. """
  1040. Compute Hessenberg form of a matrix.
  1041. The Hessenberg decomposition is::
  1042. A = Q H Q^H
  1043. where `Q` is unitary/orthogonal and `H` has only zero elements below
  1044. the first sub-diagonal.
  1045. Parameters
  1046. ----------
  1047. a : (M, M) array_like
  1048. Matrix to bring into Hessenberg form.
  1049. calc_q : bool, optional
  1050. Whether to compute the transformation matrix. Default is False.
  1051. overwrite_a : bool, optional
  1052. Whether to overwrite `a`; may improve performance.
  1053. Default is False.
  1054. check_finite : bool, optional
  1055. Whether to check that the input matrix contains only finite numbers.
  1056. Disabling may give a performance gain, but may result in problems
  1057. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  1058. Returns
  1059. -------
  1060. H : (M, M) ndarray
  1061. Hessenberg form of `a`.
  1062. Q : (M, M) ndarray
  1063. Unitary/orthogonal similarity transformation matrix ``A = Q H Q^H``.
  1064. Only returned if ``calc_q=True``.
  1065. Examples
  1066. --------
  1067. >>> from scipy.linalg import hessenberg
  1068. >>> A = np.array([[2, 5, 8, 7], [5, 2, 2, 8], [7, 5, 6, 6], [5, 4, 4, 8]])
  1069. >>> H, Q = hessenberg(A, calc_q=True)
  1070. >>> H
  1071. array([[ 2. , -11.65843866, 1.42005301, 0.25349066],
  1072. [ -9.94987437, 14.53535354, -5.31022304, 2.43081618],
  1073. [ 0. , -1.83299243, 0.38969961, -0.51527034],
  1074. [ 0. , 0. , -3.83189513, 1.07494686]])
  1075. >>> np.allclose(Q @ H @ Q.conj().T - A, np.zeros((4, 4)))
  1076. True
  1077. """
  1078. a1 = _asarray_validated(a, check_finite=check_finite)
  1079. if len(a1.shape) != 2 or (a1.shape[0] != a1.shape[1]):
  1080. raise ValueError('expected square matrix')
  1081. overwrite_a = overwrite_a or (_datacopied(a1, a))
  1082. # if 2x2 or smaller: already in Hessenberg
  1083. if a1.shape[0] <= 2:
  1084. if calc_q:
  1085. return a1, numpy.eye(a1.shape[0])
  1086. return a1
  1087. gehrd, gebal, gehrd_lwork = get_lapack_funcs(('gehrd', 'gebal',
  1088. 'gehrd_lwork'), (a1,))
  1089. ba, lo, hi, pivscale, info = gebal(a1, permute=0, overwrite_a=overwrite_a)
  1090. _check_info(info, 'gebal (hessenberg)', positive=False)
  1091. n = len(a1)
  1092. lwork = _compute_lwork(gehrd_lwork, ba.shape[0], lo=lo, hi=hi)
  1093. hq, tau, info = gehrd(ba, lo=lo, hi=hi, lwork=lwork, overwrite_a=1)
  1094. _check_info(info, 'gehrd (hessenberg)', positive=False)
  1095. h = numpy.triu(hq, -1)
  1096. if not calc_q:
  1097. return h
  1098. # use orghr/unghr to compute q
  1099. orghr, orghr_lwork = get_lapack_funcs(('orghr', 'orghr_lwork'), (a1,))
  1100. lwork = _compute_lwork(orghr_lwork, n, lo=lo, hi=hi)
  1101. q, info = orghr(a=hq, tau=tau, lo=lo, hi=hi, lwork=lwork, overwrite_a=1)
  1102. _check_info(info, 'orghr (hessenberg)', positive=False)
  1103. return h, q
  1104. def cdf2rdf(w, v):
  1105. """
  1106. Converts complex eigenvalues ``w`` and eigenvectors ``v`` to real
  1107. eigenvalues in a block diagonal form ``wr`` and the associated real
  1108. eigenvectors ``vr``, such that::
  1109. vr @ wr = X @ vr
  1110. continues to hold, where ``X`` is the original array for which ``w`` and
  1111. ``v`` are the eigenvalues and eigenvectors.
  1112. .. versionadded:: 1.1.0
  1113. Parameters
  1114. ----------
  1115. w : (..., M) array_like
  1116. Complex or real eigenvalues, an array or stack of arrays
  1117. Conjugate pairs must not be interleaved, else the wrong result
  1118. will be produced. So ``[1+1j, 1, 1-1j]`` will give a correct result, but
  1119. ``[1+1j, 2+1j, 1-1j, 2-1j]`` will not.
  1120. v : (..., M, M) array_like
  1121. Complex or real eigenvectors, a square array or stack of square arrays.
  1122. Returns
  1123. -------
  1124. wr : (..., M, M) ndarray
  1125. Real diagonal block form of eigenvalues
  1126. vr : (..., M, M) ndarray
  1127. Real eigenvectors associated with ``wr``
  1128. See Also
  1129. --------
  1130. eig : Eigenvalues and right eigenvectors for non-symmetric arrays
  1131. rsf2csf : Convert real Schur form to complex Schur form
  1132. Notes
  1133. -----
  1134. ``w``, ``v`` must be the eigenstructure for some *real* matrix ``X``.
  1135. For example, obtained by ``w, v = scipy.linalg.eig(X)`` or
  1136. ``w, v = numpy.linalg.eig(X)`` in which case ``X`` can also represent
  1137. stacked arrays.
  1138. .. versionadded:: 1.1.0
  1139. Examples
  1140. --------
  1141. >>> X = np.array([[1, 2, 3], [0, 4, 5], [0, -5, 4]])
  1142. >>> X
  1143. array([[ 1, 2, 3],
  1144. [ 0, 4, 5],
  1145. [ 0, -5, 4]])
  1146. >>> from scipy import linalg
  1147. >>> w, v = linalg.eig(X)
  1148. >>> w
  1149. array([ 1.+0.j, 4.+5.j, 4.-5.j])
  1150. >>> v
  1151. array([[ 1.00000+0.j , -0.01906-0.40016j, -0.01906+0.40016j],
  1152. [ 0.00000+0.j , 0.00000-0.64788j, 0.00000+0.64788j],
  1153. [ 0.00000+0.j , 0.64788+0.j , 0.64788-0.j ]])
  1154. >>> wr, vr = linalg.cdf2rdf(w, v)
  1155. >>> wr
  1156. array([[ 1., 0., 0.],
  1157. [ 0., 4., 5.],
  1158. [ 0., -5., 4.]])
  1159. >>> vr
  1160. array([[ 1. , 0.40016, -0.01906],
  1161. [ 0. , 0.64788, 0. ],
  1162. [ 0. , 0. , 0.64788]])
  1163. >>> vr @ wr
  1164. array([[ 1. , 1.69593, 1.9246 ],
  1165. [ 0. , 2.59153, 3.23942],
  1166. [ 0. , -3.23942, 2.59153]])
  1167. >>> X @ vr
  1168. array([[ 1. , 1.69593, 1.9246 ],
  1169. [ 0. , 2.59153, 3.23942],
  1170. [ 0. , -3.23942, 2.59153]])
  1171. """
  1172. w, v = _asarray_validated(w), _asarray_validated(v)
  1173. # check dimensions
  1174. if w.ndim < 1:
  1175. raise ValueError('expected w to be at least one-dimensional')
  1176. if v.ndim < 2:
  1177. raise ValueError('expected v to be at least two-dimensional')
  1178. if v.ndim != w.ndim + 1:
  1179. raise ValueError('expected eigenvectors array to have exactly one '
  1180. 'dimension more than eigenvalues array')
  1181. # check shapes
  1182. n = w.shape[-1]
  1183. M = w.shape[:-1]
  1184. if v.shape[-2] != v.shape[-1]:
  1185. raise ValueError('expected v to be a square matrix or stacked square '
  1186. 'matrices: v.shape[-2] = v.shape[-1]')
  1187. if v.shape[-1] != n:
  1188. raise ValueError('expected the same number of eigenvalues as '
  1189. 'eigenvectors')
  1190. # get indices for each first pair of complex eigenvalues
  1191. complex_mask = iscomplex(w)
  1192. n_complex = complex_mask.sum(axis=-1)
  1193. # check if all complex eigenvalues have conjugate pairs
  1194. if not (n_complex % 2 == 0).all():
  1195. raise ValueError('expected complex-conjugate pairs of eigenvalues')
  1196. # find complex indices
  1197. idx = nonzero(complex_mask)
  1198. idx_stack = idx[:-1]
  1199. idx_elem = idx[-1]
  1200. # filter them to conjugate indices, assuming pairs are not interleaved
  1201. j = idx_elem[0::2]
  1202. k = idx_elem[1::2]
  1203. stack_ind = ()
  1204. for i in idx_stack:
  1205. # should never happen, assuming nonzero orders by the last axis
  1206. assert (i[0::2] == i[1::2]).all(), "Conjugate pair spanned different arrays!"
  1207. stack_ind += (i[0::2],)
  1208. # all eigenvalues to diagonal form
  1209. wr = zeros(M + (n, n), dtype=w.real.dtype)
  1210. di = range(n)
  1211. wr[..., di, di] = w.real
  1212. # complex eigenvalues to real block diagonal form
  1213. wr[stack_ind + (j, k)] = w[stack_ind + (j,)].imag
  1214. wr[stack_ind + (k, j)] = w[stack_ind + (k,)].imag
  1215. # compute real eigenvectors associated with real block diagonal eigenvalues
  1216. u = zeros(M + (n, n), dtype=numpy.cdouble)
  1217. u[..., di, di] = 1.0
  1218. u[stack_ind + (j, j)] = 0.5j
  1219. u[stack_ind + (j, k)] = 0.5
  1220. u[stack_ind + (k, j)] = -0.5j
  1221. u[stack_ind + (k, k)] = 0.5
  1222. # multipy matrices v and u (equivalent to v @ u)
  1223. vr = einsum('...ij,...jk->...ik', v, u).real
  1224. return wr, vr