eigen.py 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877
  1. #!/usr/bin/python
  2. # -*- coding: utf-8 -*-
  3. ##################################################################################################
  4. # module for the eigenvalue problem
  5. # Copyright 2013 Timo Hartmann (thartmann15 at gmail.com)
  6. #
  7. # todo:
  8. # - implement balancing
  9. # - agressive early deflation
  10. #
  11. ##################################################################################################
  12. """
  13. The eigenvalue problem
  14. ----------------------
  15. This file contains routines for the eigenvalue problem.
  16. high level routines:
  17. hessenberg : reduction of a real or complex square matrix to upper Hessenberg form
  18. schur : reduction of a real or complex square matrix to upper Schur form
  19. eig : eigenvalues and eigenvectors of a real or complex square matrix
  20. low level routines:
  21. hessenberg_reduce_0 : reduction of a real or complex square matrix to upper Hessenberg form
  22. hessenberg_reduce_1 : auxiliary routine to hessenberg_reduce_0
  23. qr_step : a single implicitly shifted QR step for an upper Hessenberg matrix
  24. hessenberg_qr : Schur decomposition of an upper Hessenberg matrix
  25. eig_tr_r : right eigenvectors of an upper triangular matrix
  26. eig_tr_l : left eigenvectors of an upper triangular matrix
  27. """
  28. from ..libmp.backend import xrange
  29. class Eigen(object):
  30. pass
  31. def defun(f):
  32. setattr(Eigen, f.__name__, f)
  33. return f
  34. def hessenberg_reduce_0(ctx, A, T):
  35. """
  36. This routine computes the (upper) Hessenberg decomposition of a square matrix A.
  37. Given A, an unitary matrix Q is calculated such that
  38. Q' A Q = H and Q' Q = Q Q' = 1
  39. where H is an upper Hessenberg matrix, meaning that it only contains zeros
  40. below the first subdiagonal. Here ' denotes the hermitian transpose (i.e.
  41. transposition and conjugation).
  42. parameters:
  43. A (input/output) On input, A contains the square matrix A of
  44. dimension (n,n). On output, A contains a compressed representation
  45. of Q and H.
  46. T (output) An array of length n containing the first elements of
  47. the Householder reflectors.
  48. """
  49. # internally we work with householder reflections from the right.
  50. # let u be a row vector (i.e. u[i]=A[i,:i]). then
  51. # Q is build up by reflectors of the type (1-v'v) where v is a suitable
  52. # modification of u. these reflectors are applyed to A from the right.
  53. # because we work with reflectors from the right we have to start with
  54. # the bottom row of A and work then upwards (this corresponds to
  55. # some kind of RQ decomposition).
  56. # the first part of the vectors v (i.e. A[i,:(i-1)]) are stored as row vectors
  57. # in the lower left part of A (excluding the diagonal and subdiagonal).
  58. # the last entry of v is stored in T.
  59. # the upper right part of A (including diagonal and subdiagonal) becomes H.
  60. n = A.rows
  61. if n <= 2: return
  62. for i in xrange(n-1, 1, -1):
  63. # scale the vector
  64. scale = 0
  65. for k in xrange(0, i):
  66. scale += abs(ctx.re(A[i,k])) + abs(ctx.im(A[i,k]))
  67. scale_inv = 0
  68. if scale != 0:
  69. scale_inv = 1 / scale
  70. if scale == 0 or ctx.isinf(scale_inv):
  71. # sadly there are floating point numbers not equal to zero whose reciprocal is infinity
  72. T[i] = 0
  73. A[i,i-1] = 0
  74. continue
  75. # calculate parameters for housholder transformation
  76. H = 0
  77. for k in xrange(0, i):
  78. A[i,k] *= scale_inv
  79. rr = ctx.re(A[i,k])
  80. ii = ctx.im(A[i,k])
  81. H += rr * rr + ii * ii
  82. F = A[i,i-1]
  83. f = abs(F)
  84. G = ctx.sqrt(H)
  85. A[i,i-1] = - G * scale
  86. if f == 0:
  87. T[i] = G
  88. else:
  89. ff = F / f
  90. T[i] = F + G * ff
  91. A[i,i-1] *= ff
  92. H += G * f
  93. H = 1 / ctx.sqrt(H)
  94. T[i] *= H
  95. for k in xrange(0, i - 1):
  96. A[i,k] *= H
  97. for j in xrange(0, i):
  98. # apply housholder transformation (from right)
  99. G = ctx.conj(T[i]) * A[j,i-1]
  100. for k in xrange(0, i-1):
  101. G += ctx.conj(A[i,k]) * A[j,k]
  102. A[j,i-1] -= G * T[i]
  103. for k in xrange(0, i-1):
  104. A[j,k] -= G * A[i,k]
  105. for j in xrange(0, n):
  106. # apply housholder transformation (from left)
  107. G = T[i] * A[i-1,j]
  108. for k in xrange(0, i-1):
  109. G += A[i,k] * A[k,j]
  110. A[i-1,j] -= G * ctx.conj(T[i])
  111. for k in xrange(0, i-1):
  112. A[k,j] -= G * ctx.conj(A[i,k])
  113. def hessenberg_reduce_1(ctx, A, T):
  114. """
  115. This routine forms the unitary matrix Q described in hessenberg_reduce_0.
  116. parameters:
  117. A (input/output) On input, A is the same matrix as delivered by
  118. hessenberg_reduce_0. On output, A is set to Q.
  119. T (input) On input, T is the same array as delivered by hessenberg_reduce_0.
  120. """
  121. n = A.rows
  122. if n == 1:
  123. A[0,0] = 1
  124. return
  125. A[0,0] = A[1,1] = 1
  126. A[0,1] = A[1,0] = 0
  127. for i in xrange(2, n):
  128. if T[i] != 0:
  129. for j in xrange(0, i):
  130. G = T[i] * A[i-1,j]
  131. for k in xrange(0, i-1):
  132. G += A[i,k] * A[k,j]
  133. A[i-1,j] -= G * ctx.conj(T[i])
  134. for k in xrange(0, i-1):
  135. A[k,j] -= G * ctx.conj(A[i,k])
  136. A[i,i] = 1
  137. for j in xrange(0, i):
  138. A[j,i] = A[i,j] = 0
  139. @defun
  140. def hessenberg(ctx, A, overwrite_a = False):
  141. """
  142. This routine computes the Hessenberg decomposition of a square matrix A.
  143. Given A, an unitary matrix Q is determined such that
  144. Q' A Q = H and Q' Q = Q Q' = 1
  145. where H is an upper right Hessenberg matrix. Here ' denotes the hermitian
  146. transpose (i.e. transposition and conjugation).
  147. input:
  148. A : a real or complex square matrix
  149. overwrite_a : if true, allows modification of A which may improve
  150. performance. if false, A is not modified.
  151. output:
  152. Q : an unitary matrix
  153. H : an upper right Hessenberg matrix
  154. example:
  155. >>> from mpmath import mp
  156. >>> A = mp.matrix([[3, -1, 2], [2, 5, -5], [-2, -3, 7]])
  157. >>> Q, H = mp.hessenberg(A)
  158. >>> mp.nprint(H, 3) # doctest:+SKIP
  159. [ 3.15 2.23 4.44]
  160. [-0.769 4.85 3.05]
  161. [ 0.0 3.61 7.0]
  162. >>> print(mp.chop(A - Q * H * Q.transpose_conj()))
  163. [0.0 0.0 0.0]
  164. [0.0 0.0 0.0]
  165. [0.0 0.0 0.0]
  166. return value: (Q, H)
  167. """
  168. n = A.rows
  169. if n == 1:
  170. return (ctx.matrix([[1]]), A)
  171. if not overwrite_a:
  172. A = A.copy()
  173. T = ctx.matrix(n, 1)
  174. hessenberg_reduce_0(ctx, A, T)
  175. Q = A.copy()
  176. hessenberg_reduce_1(ctx, Q, T)
  177. for x in xrange(n):
  178. for y in xrange(x+2, n):
  179. A[y,x] = 0
  180. return Q, A
  181. ###########################################################################
  182. def qr_step(ctx, n0, n1, A, Q, shift):
  183. """
  184. This subroutine executes a single implicitly shifted QR step applied to an
  185. upper Hessenberg matrix A. Given A and shift as input, first an QR
  186. decomposition is calculated:
  187. Q R = A - shift * 1 .
  188. The output is then following matrix:
  189. R Q + shift * 1
  190. parameters:
  191. n0, n1 (input) Two integers which specify the submatrix A[n0:n1,n0:n1]
  192. on which this subroutine operators. The subdiagonal elements
  193. to the left and below this submatrix must be deflated (i.e. zero).
  194. following restriction is imposed: n1>=n0+2
  195. A (input/output) On input, A is an upper Hessenberg matrix.
  196. On output, A is replaced by "R Q + shift * 1"
  197. Q (input/output) The parameter Q is multiplied by the unitary matrix
  198. Q arising from the QR decomposition. Q can also be false, in which
  199. case the unitary matrix Q is not computated.
  200. shift (input) a complex number specifying the shift. idealy close to an
  201. eigenvalue of the bottemmost part of the submatrix A[n0:n1,n0:n1].
  202. references:
  203. Stoer, Bulirsch - Introduction to Numerical Analysis.
  204. Kresser : Numerical Methods for General and Structured Eigenvalue Problems
  205. """
  206. # implicitly shifted and bulge chasing is explained at p.398/399 in "Stoer, Bulirsch - Introduction to Numerical Analysis"
  207. # for bulge chasing see also "Watkins - The Matrix Eigenvalue Problem" sec.4.5,p.173
  208. # the Givens rotation we used is determined as follows: let c,s be two complex
  209. # numbers. then we have following relation:
  210. #
  211. # v = sqrt(|c|^2 + |s|^2)
  212. #
  213. # 1/v [ c~ s~] [c] = [v]
  214. # [-s c ] [s] [0]
  215. #
  216. # the matrix on the left is our Givens rotation.
  217. n = A.rows
  218. # first step
  219. # calculate givens rotation
  220. c = A[n0 ,n0] - shift
  221. s = A[n0+1,n0]
  222. v = ctx.hypot(ctx.hypot(ctx.re(c), ctx.im(c)), ctx.hypot(ctx.re(s), ctx.im(s)))
  223. if v == 0:
  224. v = 1
  225. c = 1
  226. s = 0
  227. else:
  228. c /= v
  229. s /= v
  230. cc = ctx.conj(c)
  231. cs = ctx.conj(s)
  232. for k in xrange(n0, n):
  233. # apply givens rotation from the left
  234. x = A[n0 ,k]
  235. y = A[n0+1,k]
  236. A[n0 ,k] = cc * x + cs * y
  237. A[n0+1,k] = c * y - s * x
  238. for k in xrange(min(n1, n0+3)):
  239. # apply givens rotation from the right
  240. x = A[k,n0 ]
  241. y = A[k,n0+1]
  242. A[k,n0 ] = c * x + s * y
  243. A[k,n0+1] = cc * y - cs * x
  244. if not isinstance(Q, bool):
  245. for k in xrange(n):
  246. # eigenvectors
  247. x = Q[k,n0 ]
  248. y = Q[k,n0+1]
  249. Q[k,n0 ] = c * x + s * y
  250. Q[k,n0+1] = cc * y - cs * x
  251. # chase the bulge
  252. for j in xrange(n0, n1 - 2):
  253. # calculate givens rotation
  254. c = A[j+1,j]
  255. s = A[j+2,j]
  256. v = ctx.hypot(ctx.hypot(ctx.re(c), ctx.im(c)), ctx.hypot(ctx.re(s), ctx.im(s)))
  257. if v == 0:
  258. A[j+1,j] = 0
  259. v = 1
  260. c = 1
  261. s = 0
  262. else:
  263. A[j+1,j] = v
  264. c /= v
  265. s /= v
  266. A[j+2,j] = 0
  267. cc = ctx.conj(c)
  268. cs = ctx.conj(s)
  269. for k in xrange(j+1, n):
  270. # apply givens rotation from the left
  271. x = A[j+1,k]
  272. y = A[j+2,k]
  273. A[j+1,k] = cc * x + cs * y
  274. A[j+2,k] = c * y - s * x
  275. for k in xrange(0, min(n1, j+4)):
  276. # apply givens rotation from the right
  277. x = A[k,j+1]
  278. y = A[k,j+2]
  279. A[k,j+1] = c * x + s * y
  280. A[k,j+2] = cc * y - cs * x
  281. if not isinstance(Q, bool):
  282. for k in xrange(0, n):
  283. # eigenvectors
  284. x = Q[k,j+1]
  285. y = Q[k,j+2]
  286. Q[k,j+1] = c * x + s * y
  287. Q[k,j+2] = cc * y - cs * x
  288. def hessenberg_qr(ctx, A, Q):
  289. """
  290. This routine computes the Schur decomposition of an upper Hessenberg matrix A.
  291. Given A, an unitary matrix Q is determined such that
  292. Q' A Q = R and Q' Q = Q Q' = 1
  293. where R is an upper right triangular matrix. Here ' denotes the hermitian
  294. transpose (i.e. transposition and conjugation).
  295. parameters:
  296. A (input/output) On input, A contains an upper Hessenberg matrix.
  297. On output, A is replace by the upper right triangluar matrix R.
  298. Q (input/output) The parameter Q is multiplied by the unitary
  299. matrix Q arising from the Schur decomposition. Q can also be
  300. false, in which case the unitary matrix Q is not computated.
  301. """
  302. n = A.rows
  303. norm = 0
  304. for x in xrange(n):
  305. for y in xrange(min(x+2, n)):
  306. norm += ctx.re(A[y,x]) ** 2 + ctx.im(A[y,x]) ** 2
  307. norm = ctx.sqrt(norm) / n
  308. if norm == 0:
  309. return
  310. n0 = 0
  311. n1 = n
  312. eps = ctx.eps / (100 * n)
  313. maxits = ctx.dps * 4
  314. its = totalits = 0
  315. while 1:
  316. # kressner p.32 algo 3
  317. # the active submatrix is A[n0:n1,n0:n1]
  318. k = n0
  319. while k + 1 < n1:
  320. s = abs(ctx.re(A[k,k])) + abs(ctx.im(A[k,k])) + abs(ctx.re(A[k+1,k+1])) + abs(ctx.im(A[k+1,k+1]))
  321. if s < eps * norm:
  322. s = norm
  323. if abs(A[k+1,k]) < eps * s:
  324. break
  325. k += 1
  326. if k + 1 < n1:
  327. # deflation found at position (k+1, k)
  328. A[k+1,k] = 0
  329. n0 = k + 1
  330. its = 0
  331. if n0 + 1 >= n1:
  332. # block of size at most two has converged
  333. n0 = 0
  334. n1 = k + 1
  335. if n1 < 2:
  336. # QR algorithm has converged
  337. return
  338. else:
  339. if (its % 30) == 10:
  340. # exceptional shift
  341. shift = A[n1-1,n1-2]
  342. elif (its % 30) == 20:
  343. # exceptional shift
  344. shift = abs(A[n1-1,n1-2])
  345. elif (its % 30) == 29:
  346. # exceptional shift
  347. shift = norm
  348. else:
  349. # A = [ a b ] det(x-A)=x*x-x*tr(A)+det(A)
  350. # [ c d ]
  351. #
  352. # eigenvalues bad: (tr(A)+sqrt((tr(A))**2-4*det(A)))/2
  353. # bad because of cancellation if |c| is small and |a-d| is small, too.
  354. #
  355. # eigenvalues good: (a+d+sqrt((a-d)**2+4*b*c))/2
  356. t = A[n1-2,n1-2] + A[n1-1,n1-1]
  357. s = (A[n1-1,n1-1] - A[n1-2,n1-2]) ** 2 + 4 * A[n1-1,n1-2] * A[n1-2,n1-1]
  358. if ctx.re(s) > 0:
  359. s = ctx.sqrt(s)
  360. else:
  361. s = ctx.sqrt(-s) * 1j
  362. a = (t + s) / 2
  363. b = (t - s) / 2
  364. if abs(A[n1-1,n1-1] - a) > abs(A[n1-1,n1-1] - b):
  365. shift = b
  366. else:
  367. shift = a
  368. its += 1
  369. totalits += 1
  370. qr_step(ctx, n0, n1, A, Q, shift)
  371. if its > maxits:
  372. raise RuntimeError("qr: failed to converge after %d steps" % its)
  373. @defun
  374. def schur(ctx, A, overwrite_a = False):
  375. """
  376. This routine computes the Schur decomposition of a square matrix A.
  377. Given A, an unitary matrix Q is determined such that
  378. Q' A Q = R and Q' Q = Q Q' = 1
  379. where R is an upper right triangular matrix. Here ' denotes the
  380. hermitian transpose (i.e. transposition and conjugation).
  381. input:
  382. A : a real or complex square matrix
  383. overwrite_a : if true, allows modification of A which may improve
  384. performance. if false, A is not modified.
  385. output:
  386. Q : an unitary matrix
  387. R : an upper right triangular matrix
  388. return value: (Q, R)
  389. example:
  390. >>> from mpmath import mp
  391. >>> A = mp.matrix([[3, -1, 2], [2, 5, -5], [-2, -3, 7]])
  392. >>> Q, R = mp.schur(A)
  393. >>> mp.nprint(R, 3) # doctest:+SKIP
  394. [2.0 0.417 -2.53]
  395. [0.0 4.0 -4.74]
  396. [0.0 0.0 9.0]
  397. >>> print(mp.chop(A - Q * R * Q.transpose_conj()))
  398. [0.0 0.0 0.0]
  399. [0.0 0.0 0.0]
  400. [0.0 0.0 0.0]
  401. warning: The Schur decomposition is not unique.
  402. """
  403. n = A.rows
  404. if n == 1:
  405. return (ctx.matrix([[1]]), A)
  406. if not overwrite_a:
  407. A = A.copy()
  408. T = ctx.matrix(n, 1)
  409. hessenberg_reduce_0(ctx, A, T)
  410. Q = A.copy()
  411. hessenberg_reduce_1(ctx, Q, T)
  412. for x in xrange(n):
  413. for y in xrange(x + 2, n):
  414. A[y,x] = 0
  415. hessenberg_qr(ctx, A, Q)
  416. return Q, A
  417. def eig_tr_r(ctx, A):
  418. """
  419. This routine calculates the right eigenvectors of an upper right triangular matrix.
  420. input:
  421. A an upper right triangular matrix
  422. output:
  423. ER a matrix whose columns form the right eigenvectors of A
  424. return value: ER
  425. """
  426. # this subroutine is inspired by the lapack routines ctrevc.f,clatrs.f
  427. n = A.rows
  428. ER = ctx.eye(n)
  429. eps = ctx.eps
  430. unfl = ctx.ldexp(ctx.one, -ctx.prec * 30)
  431. # since mpmath effectively has no limits on the exponent, we simply scale doubles up
  432. # original double has prec*20
  433. smlnum = unfl * (n / eps)
  434. simin = 1 / ctx.sqrt(eps)
  435. rmax = 1
  436. for i in xrange(1, n):
  437. s = A[i,i]
  438. smin = max(eps * abs(s), smlnum)
  439. for j in xrange(i - 1, -1, -1):
  440. r = 0
  441. for k in xrange(j + 1, i + 1):
  442. r += A[j,k] * ER[k,i]
  443. t = A[j,j] - s
  444. if abs(t) < smin:
  445. t = smin
  446. r = -r / t
  447. ER[j,i] = r
  448. rmax = max(rmax, abs(r))
  449. if rmax > simin:
  450. for k in xrange(j, i+1):
  451. ER[k,i] /= rmax
  452. rmax = 1
  453. if rmax != 1:
  454. for k in xrange(0, i + 1):
  455. ER[k,i] /= rmax
  456. return ER
  457. def eig_tr_l(ctx, A):
  458. """
  459. This routine calculates the left eigenvectors of an upper right triangular matrix.
  460. input:
  461. A an upper right triangular matrix
  462. output:
  463. EL a matrix whose rows form the left eigenvectors of A
  464. return value: EL
  465. """
  466. n = A.rows
  467. EL = ctx.eye(n)
  468. eps = ctx.eps
  469. unfl = ctx.ldexp(ctx.one, -ctx.prec * 30)
  470. # since mpmath effectively has no limits on the exponent, we simply scale doubles up
  471. # original double has prec*20
  472. smlnum = unfl * (n / eps)
  473. simin = 1 / ctx.sqrt(eps)
  474. rmax = 1
  475. for i in xrange(0, n - 1):
  476. s = A[i,i]
  477. smin = max(eps * abs(s), smlnum)
  478. for j in xrange(i + 1, n):
  479. r = 0
  480. for k in xrange(i, j):
  481. r += EL[i,k] * A[k,j]
  482. t = A[j,j] - s
  483. if abs(t) < smin:
  484. t = smin
  485. r = -r / t
  486. EL[i,j] = r
  487. rmax = max(rmax, abs(r))
  488. if rmax > simin:
  489. for k in xrange(i, j + 1):
  490. EL[i,k] /= rmax
  491. rmax = 1
  492. if rmax != 1:
  493. for k in xrange(i, n):
  494. EL[i,k] /= rmax
  495. return EL
  496. @defun
  497. def eig(ctx, A, left = False, right = True, overwrite_a = False):
  498. """
  499. This routine computes the eigenvalues and optionally the left and right
  500. eigenvectors of a square matrix A. Given A, a vector E and matrices ER
  501. and EL are calculated such that
  502. A ER[:,i] = E[i] ER[:,i]
  503. EL[i,:] A = EL[i,:] E[i]
  504. E contains the eigenvalues of A. The columns of ER contain the right eigenvectors
  505. of A whereas the rows of EL contain the left eigenvectors.
  506. input:
  507. A : a real or complex square matrix of shape (n, n)
  508. left : if true, the left eigenvectors are calculated.
  509. right : if true, the right eigenvectors are calculated.
  510. overwrite_a : if true, allows modification of A which may improve
  511. performance. if false, A is not modified.
  512. output:
  513. E : a list of length n containing the eigenvalues of A.
  514. ER : a matrix whose columns contain the right eigenvectors of A.
  515. EL : a matrix whose rows contain the left eigenvectors of A.
  516. return values:
  517. E if left and right are both false.
  518. (E, ER) if right is true and left is false.
  519. (E, EL) if left is true and right is false.
  520. (E, EL, ER) if left and right are true.
  521. examples:
  522. >>> from mpmath import mp
  523. >>> A = mp.matrix([[3, -1, 2], [2, 5, -5], [-2, -3, 7]])
  524. >>> E, ER = mp.eig(A)
  525. >>> print(mp.chop(A * ER[:,0] - E[0] * ER[:,0]))
  526. [0.0]
  527. [0.0]
  528. [0.0]
  529. >>> E, EL, ER = mp.eig(A,left = True, right = True)
  530. >>> E, EL, ER = mp.eig_sort(E, EL, ER)
  531. >>> mp.nprint(E)
  532. [2.0, 4.0, 9.0]
  533. >>> print(mp.chop(A * ER[:,0] - E[0] * ER[:,0]))
  534. [0.0]
  535. [0.0]
  536. [0.0]
  537. >>> print(mp.chop( EL[0,:] * A - EL[0,:] * E[0]))
  538. [0.0 0.0 0.0]
  539. warning:
  540. - If there are multiple eigenvalues, the eigenvectors do not necessarily
  541. span the whole vectorspace, i.e. ER and EL may have not full rank.
  542. Furthermore in that case the eigenvectors are numerical ill-conditioned.
  543. - In the general case the eigenvalues have no natural order.
  544. see also:
  545. - eigh (or eigsy, eighe) for the symmetric eigenvalue problem.
  546. - eig_sort for sorting of eigenvalues and eigenvectors
  547. """
  548. n = A.rows
  549. if n == 1:
  550. if left and (not right):
  551. return ([A[0]], ctx.matrix([[1]]))
  552. if right and (not left):
  553. return ([A[0]], ctx.matrix([[1]]))
  554. return ([A[0]], ctx.matrix([[1]]), ctx.matrix([[1]]))
  555. if not overwrite_a:
  556. A = A.copy()
  557. T = ctx.zeros(n, 1)
  558. hessenberg_reduce_0(ctx, A, T)
  559. if left or right:
  560. Q = A.copy()
  561. hessenberg_reduce_1(ctx, Q, T)
  562. else:
  563. Q = False
  564. for x in xrange(n):
  565. for y in xrange(x + 2, n):
  566. A[y,x] = 0
  567. hessenberg_qr(ctx, A, Q)
  568. E = [0 for i in xrange(n)]
  569. for i in xrange(n):
  570. E[i] = A[i,i]
  571. if not (left or right):
  572. return E
  573. if left:
  574. EL = eig_tr_l(ctx, A)
  575. EL = EL * Q.transpose_conj()
  576. if right:
  577. ER = eig_tr_r(ctx, A)
  578. ER = Q * ER
  579. if left and (not right):
  580. return (E, EL)
  581. if right and (not left):
  582. return (E, ER)
  583. return (E, EL, ER)
  584. @defun
  585. def eig_sort(ctx, E, EL = False, ER = False, f = "real"):
  586. """
  587. This routine sorts the eigenvalues and eigenvectors delivered by ``eig``.
  588. parameters:
  589. E : the eigenvalues as delivered by eig
  590. EL : the left eigenvectors as delivered by eig, or false
  591. ER : the right eigenvectors as delivered by eig, or false
  592. f : either a string ("real" sort by increasing real part, "imag" sort by
  593. increasing imag part, "abs" sort by absolute value) or a function
  594. mapping complexs to the reals, i.e. ``f = lambda x: -mp.re(x) ``
  595. would sort the eigenvalues by decreasing real part.
  596. return values:
  597. E if EL and ER are both false.
  598. (E, ER) if ER is not false and left is false.
  599. (E, EL) if EL is not false and right is false.
  600. (E, EL, ER) if EL and ER are not false.
  601. example:
  602. >>> from mpmath import mp
  603. >>> A = mp.matrix([[3, -1, 2], [2, 5, -5], [-2, -3, 7]])
  604. >>> E, EL, ER = mp.eig(A,left = True, right = True)
  605. >>> E, EL, ER = mp.eig_sort(E, EL, ER)
  606. >>> mp.nprint(E)
  607. [2.0, 4.0, 9.0]
  608. >>> E, EL, ER = mp.eig_sort(E, EL, ER,f = lambda x: -mp.re(x))
  609. >>> mp.nprint(E)
  610. [9.0, 4.0, 2.0]
  611. >>> print(mp.chop(A * ER[:,0] - E[0] * ER[:,0]))
  612. [0.0]
  613. [0.0]
  614. [0.0]
  615. >>> print(mp.chop( EL[0,:] * A - EL[0,:] * E[0]))
  616. [0.0 0.0 0.0]
  617. """
  618. if isinstance(f, str):
  619. if f == "real":
  620. f = ctx.re
  621. elif f == "imag":
  622. f = ctx.im
  623. elif f == "abs":
  624. f = abs
  625. else:
  626. raise RuntimeError("unknown function %s" % f)
  627. n = len(E)
  628. # Sort eigenvalues (bubble-sort)
  629. for i in xrange(n):
  630. imax = i
  631. s = f(E[i]) # s is the current maximal element
  632. for j in xrange(i + 1, n):
  633. c = f(E[j])
  634. if c < s:
  635. s = c
  636. imax = j
  637. if imax != i:
  638. # swap eigenvalues
  639. z = E[i]
  640. E[i] = E[imax]
  641. E[imax] = z
  642. if not isinstance(EL, bool):
  643. for j in xrange(n):
  644. z = EL[i,j]
  645. EL[i,j] = EL[imax,j]
  646. EL[imax,j] = z
  647. if not isinstance(ER, bool):
  648. for j in xrange(n):
  649. z = ER[j,i]
  650. ER[j,i] = ER[j,imax]
  651. ER[j,imax] = z
  652. if isinstance(EL, bool) and isinstance(ER, bool):
  653. return E
  654. if isinstance(EL, bool) and not(isinstance(ER, bool)):
  655. return (E, ER)
  656. if isinstance(ER, bool) and not(isinstance(EL, bool)):
  657. return (E, EL)
  658. return (E, EL, ER)