eigen_symmetric.py 57 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807
  1. #!/usr/bin/python
  2. # -*- coding: utf-8 -*-
  3. ##################################################################################################
  4. # module for the symmetric eigenvalue problem
  5. # Copyright 2013 Timo Hartmann (thartmann15 at gmail.com)
  6. #
  7. # todo:
  8. # - implement balancing
  9. #
  10. ##################################################################################################
  11. """
  12. The symmetric eigenvalue problem.
  13. ---------------------------------
  14. This file contains routines for the symmetric eigenvalue problem.
  15. high level routines:
  16. eigsy : real symmetric (ordinary) eigenvalue problem
  17. eighe : complex hermitian (ordinary) eigenvalue problem
  18. eigh : unified interface for eigsy and eighe
  19. svd_r : singular value decomposition for real matrices
  20. svd_c : singular value decomposition for complex matrices
  21. svd : unified interface for svd_r and svd_c
  22. low level routines:
  23. r_sy_tridiag : reduction of real symmetric matrix to real symmetric tridiagonal matrix
  24. c_he_tridiag_0 : reduction of complex hermitian matrix to real symmetric tridiagonal matrix
  25. c_he_tridiag_1 : auxiliary routine to c_he_tridiag_0
  26. c_he_tridiag_2 : auxiliary routine to c_he_tridiag_0
  27. tridiag_eigen : solves the real symmetric tridiagonal matrix eigenvalue problem
  28. svd_r_raw : raw singular value decomposition for real matrices
  29. svd_c_raw : raw singular value decomposition for complex matrices
  30. """
  31. from ..libmp.backend import xrange
  32. from .eigen import defun
  33. def r_sy_tridiag(ctx, A, D, E, calc_ev = True):
  34. """
  35. This routine transforms a real symmetric matrix A to a real symmetric
  36. tridiagonal matrix T using an orthogonal similarity transformation:
  37. Q' * A * Q = T (here ' denotes the matrix transpose).
  38. The orthogonal matrix Q is build up from Householder reflectors.
  39. parameters:
  40. A (input/output) On input, A contains the real symmetric matrix of
  41. dimension (n,n). On output, if calc_ev is true, A contains the
  42. orthogonal matrix Q, otherwise A is destroyed.
  43. D (output) real array of length n, contains the diagonal elements
  44. of the tridiagonal matrix
  45. E (output) real array of length n, contains the offdiagonal elements
  46. of the tridiagonal matrix in E[0:(n-1)] where is the dimension of
  47. the matrix A. E[n-1] is undefined.
  48. calc_ev (input) If calc_ev is true, this routine explicitly calculates the
  49. orthogonal matrix Q which is then returned in A. If calc_ev is
  50. false, Q is not explicitly calculated resulting in a shorter run time.
  51. This routine is a python translation of the fortran routine tred2.f in the
  52. software library EISPACK (see netlib.org) which itself is based on the algol
  53. procedure tred2 described in:
  54. - Num. Math. 11, p.181-195 (1968) by Martin, Reinsch and Wilkonson
  55. - Handbook for auto. comp., Vol II, Linear Algebra, p.212-226 (1971)
  56. For a good introduction to Householder reflections, see also
  57. Stoer, Bulirsch - Introduction to Numerical Analysis.
  58. """
  59. # note : the vector v of the i-th houshoulder reflector is stored in a[(i+1):,i]
  60. # whereas v/<v,v> is stored in a[i,(i+1):]
  61. n = A.rows
  62. for i in xrange(n - 1, 0, -1):
  63. # scale the vector
  64. scale = 0
  65. for k in xrange(0, i):
  66. scale += abs(A[k,i])
  67. scale_inv = 0
  68. if scale != 0:
  69. scale_inv = 1/scale
  70. # sadly there are floating point numbers not equal to zero whose reciprocal is infinity
  71. if i == 1 or scale == 0 or ctx.isinf(scale_inv):
  72. E[i] = A[i-1,i] # nothing to do
  73. D[i] = 0
  74. continue
  75. # calculate parameters for housholder transformation
  76. H = 0
  77. for k in xrange(0, i):
  78. A[k,i] *= scale_inv
  79. H += A[k,i] * A[k,i]
  80. F = A[i-1,i]
  81. G = ctx.sqrt(H)
  82. if F > 0:
  83. G = -G
  84. E[i] = scale * G
  85. H -= F * G
  86. A[i-1,i] = F - G
  87. F = 0
  88. # apply housholder transformation
  89. for j in xrange(0, i):
  90. if calc_ev:
  91. A[i,j] = A[j,i] / H
  92. G = 0 # calculate A*U
  93. for k in xrange(0, j + 1):
  94. G += A[k,j] * A[k,i]
  95. for k in xrange(j + 1, i):
  96. G += A[j,k] * A[k,i]
  97. E[j] = G / H # calculate P
  98. F += E[j] * A[j,i]
  99. HH = F / (2 * H)
  100. for j in xrange(0, i): # calculate reduced A
  101. F = A[j,i]
  102. G = E[j] - HH * F # calculate Q
  103. E[j] = G
  104. for k in xrange(0, j + 1):
  105. A[k,j] -= F * E[k] + G * A[k,i]
  106. D[i] = H
  107. for i in xrange(1, n): # better for compatibility
  108. E[i-1] = E[i]
  109. E[n-1] = 0
  110. if calc_ev:
  111. D[0] = 0
  112. for i in xrange(0, n):
  113. if D[i] != 0:
  114. for j in xrange(0, i): # accumulate transformation matrices
  115. G = 0
  116. for k in xrange(0, i):
  117. G += A[i,k] * A[k,j]
  118. for k in xrange(0, i):
  119. A[k,j] -= G * A[k,i]
  120. D[i] = A[i,i]
  121. A[i,i] = 1
  122. for j in xrange(0, i):
  123. A[j,i] = A[i,j] = 0
  124. else:
  125. for i in xrange(0, n):
  126. D[i] = A[i,i]
  127. def c_he_tridiag_0(ctx, A, D, E, T):
  128. """
  129. This routine transforms a complex hermitian matrix A to a real symmetric
  130. tridiagonal matrix T using an unitary similarity transformation:
  131. Q' * A * Q = T (here ' denotes the hermitian matrix transpose,
  132. i.e. transposition und conjugation).
  133. The unitary matrix Q is build up from Householder reflectors and
  134. an unitary diagonal matrix.
  135. parameters:
  136. A (input/output) On input, A contains the complex hermitian matrix
  137. of dimension (n,n). On output, A contains the unitary matrix Q
  138. in compressed form.
  139. D (output) real array of length n, contains the diagonal elements
  140. of the tridiagonal matrix.
  141. E (output) real array of length n, contains the offdiagonal elements
  142. of the tridiagonal matrix in E[0:(n-1)] where is the dimension of
  143. the matrix A. E[n-1] is undefined.
  144. T (output) complex array of length n, contains a unitary diagonal
  145. matrix.
  146. This routine is a python translation (in slightly modified form) of the fortran
  147. routine htridi.f in the software library EISPACK (see netlib.org) which itself
  148. is a complex version of the algol procedure tred1 described in:
  149. - Num. Math. 11, p.181-195 (1968) by Martin, Reinsch and Wilkonson
  150. - Handbook for auto. comp., Vol II, Linear Algebra, p.212-226 (1971)
  151. For a good introduction to Householder reflections, see also
  152. Stoer, Bulirsch - Introduction to Numerical Analysis.
  153. """
  154. n = A.rows
  155. T[n-1] = 1
  156. for i in xrange(n - 1, 0, -1):
  157. # scale the vector
  158. scale = 0
  159. for k in xrange(0, i):
  160. scale += abs(ctx.re(A[k,i])) + abs(ctx.im(A[k,i]))
  161. scale_inv = 0
  162. if scale != 0:
  163. scale_inv = 1 / scale
  164. # sadly there are floating point numbers not equal to zero whose reciprocal is infinity
  165. if scale == 0 or ctx.isinf(scale_inv):
  166. E[i] = 0
  167. D[i] = 0
  168. T[i-1] = 1
  169. continue
  170. if i == 1:
  171. F = A[i-1,i]
  172. f = abs(F)
  173. E[i] = f
  174. D[i] = 0
  175. if f != 0:
  176. T[i-1] = T[i] * F / f
  177. else:
  178. T[i-1] = T[i]
  179. continue
  180. # calculate parameters for housholder transformation
  181. H = 0
  182. for k in xrange(0, i):
  183. A[k,i] *= scale_inv
  184. rr = ctx.re(A[k,i])
  185. ii = ctx.im(A[k,i])
  186. H += rr * rr + ii * ii
  187. F = A[i-1,i]
  188. f = abs(F)
  189. G = ctx.sqrt(H)
  190. H += G * f
  191. E[i] = scale * G
  192. if f != 0:
  193. F = F / f
  194. TZ = - T[i] * F # T[i-1]=-T[i]*F, but we need T[i-1] as temporary storage
  195. G *= F
  196. else:
  197. TZ = -T[i] # T[i-1]=-T[i]
  198. A[i-1,i] += G
  199. F = 0
  200. # apply housholder transformation
  201. for j in xrange(0, i):
  202. A[i,j] = A[j,i] / H
  203. G = 0 # calculate A*U
  204. for k in xrange(0, j + 1):
  205. G += ctx.conj(A[k,j]) * A[k,i]
  206. for k in xrange(j + 1, i):
  207. G += A[j,k] * A[k,i]
  208. T[j] = G / H # calculate P
  209. F += ctx.conj(T[j]) * A[j,i]
  210. HH = F / (2 * H)
  211. for j in xrange(0, i): # calculate reduced A
  212. F = A[j,i]
  213. G = T[j] - HH * F # calculate Q
  214. T[j] = G
  215. for k in xrange(0, j + 1):
  216. A[k,j] -= ctx.conj(F) * T[k] + ctx.conj(G) * A[k,i]
  217. # as we use the lower left part for storage
  218. # we have to use the transpose of the normal formula
  219. T[i-1] = TZ
  220. D[i] = H
  221. for i in xrange(1, n): # better for compatibility
  222. E[i-1] = E[i]
  223. E[n-1] = 0
  224. D[0] = 0
  225. for i in xrange(0, n):
  226. zw = D[i]
  227. D[i] = ctx.re(A[i,i])
  228. A[i,i] = zw
  229. def c_he_tridiag_1(ctx, A, T):
  230. """
  231. This routine forms the unitary matrix Q described in c_he_tridiag_0.
  232. parameters:
  233. A (input/output) On input, A is the same matrix as delivered by
  234. c_he_tridiag_0. On output, A is set to Q.
  235. T (input) On input, T is the same array as delivered by c_he_tridiag_0.
  236. """
  237. n = A.rows
  238. for i in xrange(0, n):
  239. if A[i,i] != 0:
  240. for j in xrange(0, i):
  241. G = 0
  242. for k in xrange(0, i):
  243. G += ctx.conj(A[i,k]) * A[k,j]
  244. for k in xrange(0, i):
  245. A[k,j] -= G * A[k,i]
  246. A[i,i] = 1
  247. for j in xrange(0, i):
  248. A[j,i] = A[i,j] = 0
  249. for i in xrange(0, n):
  250. for k in xrange(0, n):
  251. A[i,k] *= T[k]
  252. def c_he_tridiag_2(ctx, A, T, B):
  253. """
  254. This routine applied the unitary matrix Q described in c_he_tridiag_0
  255. onto the the matrix B, i.e. it forms Q*B.
  256. parameters:
  257. A (input) On input, A is the same matrix as delivered by c_he_tridiag_0.
  258. T (input) On input, T is the same array as delivered by c_he_tridiag_0.
  259. B (input/output) On input, B is a complex matrix. On output B is replaced
  260. by Q*B.
  261. This routine is a python translation of the fortran routine htribk.f in the
  262. software library EISPACK (see netlib.org). See c_he_tridiag_0 for more
  263. references.
  264. """
  265. n = A.rows
  266. for i in xrange(0, n):
  267. for k in xrange(0, n):
  268. B[k,i] *= T[k]
  269. for i in xrange(0, n):
  270. if A[i,i] != 0:
  271. for j in xrange(0, n):
  272. G = 0
  273. for k in xrange(0, i):
  274. G += ctx.conj(A[i,k]) * B[k,j]
  275. for k in xrange(0, i):
  276. B[k,j] -= G * A[k,i]
  277. def tridiag_eigen(ctx, d, e, z = False):
  278. """
  279. This subroutine find the eigenvalues and the first components of the
  280. eigenvectors of a real symmetric tridiagonal matrix using the implicit
  281. QL method.
  282. parameters:
  283. d (input/output) real array of length n. on input, d contains the diagonal
  284. elements of the input matrix. on output, d contains the eigenvalues in
  285. ascending order.
  286. e (input) real array of length n. on input, e contains the offdiagonal
  287. elements of the input matrix in e[0:(n-1)]. On output, e has been
  288. destroyed.
  289. z (input/output) If z is equal to False, no eigenvectors will be computed.
  290. Otherwise on input z should have the format z[0:m,0:n] (i.e. a real or
  291. complex matrix of dimension (m,n) ). On output this matrix will be
  292. multiplied by the matrix of the eigenvectors (i.e. the columns of this
  293. matrix are the eigenvectors): z --> z*EV
  294. That means if z[i,j]={1 if j==j; 0 otherwise} on input, then on output
  295. z will contain the first m components of the eigenvectors. That means
  296. if m is equal to n, the i-th eigenvector will be z[:,i].
  297. This routine is a python translation (in slightly modified form) of the
  298. fortran routine imtql2.f in the software library EISPACK (see netlib.org)
  299. which itself is based on the algol procudure imtql2 desribed in:
  300. - num. math. 12, p. 377-383(1968) by matrin and wilkinson
  301. - modified in num. math. 15, p. 450(1970) by dubrulle
  302. - handbook for auto. comp., vol. II-linear algebra, p. 241-248 (1971)
  303. See also the routine gaussq.f in netlog.org or acm algorithm 726.
  304. """
  305. n = len(d)
  306. e[n-1] = 0
  307. iterlim = 2 * ctx.dps
  308. for l in xrange(n):
  309. j = 0
  310. while 1:
  311. m = l
  312. while 1:
  313. # look for a small subdiagonal element
  314. if m + 1 == n:
  315. break
  316. if abs(e[m]) <= ctx.eps * (abs(d[m]) + abs(d[m + 1])):
  317. break
  318. m = m + 1
  319. if m == l:
  320. break
  321. if j >= iterlim:
  322. raise RuntimeError("tridiag_eigen: no convergence to an eigenvalue after %d iterations" % iterlim)
  323. j += 1
  324. # form shift
  325. p = d[l]
  326. g = (d[l + 1] - p) / (2 * e[l])
  327. r = ctx.hypot(g, 1)
  328. if g < 0:
  329. s = g - r
  330. else:
  331. s = g + r
  332. g = d[m] - p + e[l] / s
  333. s, c, p = 1, 1, 0
  334. for i in xrange(m - 1, l - 1, -1):
  335. f = s * e[i]
  336. b = c * e[i]
  337. if abs(f) > abs(g): # this here is a slight improvement also used in gaussq.f or acm algorithm 726.
  338. c = g / f
  339. r = ctx.hypot(c, 1)
  340. e[i + 1] = f * r
  341. s = 1 / r
  342. c = c * s
  343. else:
  344. s = f / g
  345. r = ctx.hypot(s, 1)
  346. e[i + 1] = g * r
  347. c = 1 / r
  348. s = s * c
  349. g = d[i + 1] - p
  350. r = (d[i] - g) * s + 2 * c * b
  351. p = s * r
  352. d[i + 1] = g + p
  353. g = c * r - b
  354. if not isinstance(z, bool):
  355. # calculate eigenvectors
  356. for w in xrange(z.rows):
  357. f = z[w,i+1]
  358. z[w,i+1] = s * z[w,i] + c * f
  359. z[w,i ] = c * z[w,i] - s * f
  360. d[l] = d[l] - p
  361. e[l] = g
  362. e[m] = 0
  363. for ii in xrange(1, n):
  364. # sort eigenvalues and eigenvectors (bubble-sort)
  365. i = ii - 1
  366. k = i
  367. p = d[i]
  368. for j in xrange(ii, n):
  369. if d[j] >= p:
  370. continue
  371. k = j
  372. p = d[k]
  373. if k == i:
  374. continue
  375. d[k] = d[i]
  376. d[i] = p
  377. if not isinstance(z, bool):
  378. for w in xrange(z.rows):
  379. p = z[w,i]
  380. z[w,i] = z[w,k]
  381. z[w,k] = p
  382. ########################################################################################
  383. @defun
  384. def eigsy(ctx, A, eigvals_only = False, overwrite_a = False):
  385. """
  386. This routine solves the (ordinary) eigenvalue problem for a real symmetric
  387. square matrix A. Given A, an orthogonal matrix Q is calculated which
  388. diagonalizes A:
  389. Q' A Q = diag(E) and Q Q' = Q' Q = 1
  390. Here diag(E) is a diagonal matrix whose diagonal is E.
  391. ' denotes the transpose.
  392. The columns of Q are the eigenvectors of A and E contains the eigenvalues:
  393. A Q[:,i] = E[i] Q[:,i]
  394. input:
  395. A: real matrix of format (n,n) which is symmetric
  396. (i.e. A=A' or A[i,j]=A[j,i])
  397. eigvals_only: if true, calculates only the eigenvalues E.
  398. if false, calculates both eigenvectors and eigenvalues.
  399. overwrite_a: if true, allows modification of A which may improve
  400. performance. if false, A is not modified.
  401. output:
  402. E: vector of format (n). contains the eigenvalues of A in ascending order.
  403. Q: orthogonal matrix of format (n,n). contains the eigenvectors
  404. of A as columns.
  405. return value:
  406. E if eigvals_only is true
  407. (E, Q) if eigvals_only is false
  408. example:
  409. >>> from mpmath import mp
  410. >>> A = mp.matrix([[3, 2], [2, 0]])
  411. >>> E = mp.eigsy(A, eigvals_only = True)
  412. >>> print(E)
  413. [-1.0]
  414. [ 4.0]
  415. >>> A = mp.matrix([[1, 2], [2, 3]])
  416. >>> E, Q = mp.eigsy(A)
  417. >>> print(mp.chop(A * Q[:,0] - E[0] * Q[:,0]))
  418. [0.0]
  419. [0.0]
  420. see also: eighe, eigh, eig
  421. """
  422. if not overwrite_a:
  423. A = A.copy()
  424. d = ctx.zeros(A.rows, 1)
  425. e = ctx.zeros(A.rows, 1)
  426. if eigvals_only:
  427. r_sy_tridiag(ctx, A, d, e, calc_ev = False)
  428. tridiag_eigen(ctx, d, e, False)
  429. return d
  430. else:
  431. r_sy_tridiag(ctx, A, d, e, calc_ev = True)
  432. tridiag_eigen(ctx, d, e, A)
  433. return (d, A)
  434. @defun
  435. def eighe(ctx, A, eigvals_only = False, overwrite_a = False):
  436. """
  437. This routine solves the (ordinary) eigenvalue problem for a complex
  438. hermitian square matrix A. Given A, an unitary matrix Q is calculated which
  439. diagonalizes A:
  440. Q' A Q = diag(E) and Q Q' = Q' Q = 1
  441. Here diag(E) a is diagonal matrix whose diagonal is E.
  442. ' denotes the hermitian transpose (i.e. ordinary transposition and
  443. complex conjugation).
  444. The columns of Q are the eigenvectors of A and E contains the eigenvalues:
  445. A Q[:,i] = E[i] Q[:,i]
  446. input:
  447. A: complex matrix of format (n,n) which is hermitian
  448. (i.e. A=A' or A[i,j]=conj(A[j,i]))
  449. eigvals_only: if true, calculates only the eigenvalues E.
  450. if false, calculates both eigenvectors and eigenvalues.
  451. overwrite_a: if true, allows modification of A which may improve
  452. performance. if false, A is not modified.
  453. output:
  454. E: vector of format (n). contains the eigenvalues of A in ascending order.
  455. Q: unitary matrix of format (n,n). contains the eigenvectors
  456. of A as columns.
  457. return value:
  458. E if eigvals_only is true
  459. (E, Q) if eigvals_only is false
  460. example:
  461. >>> from mpmath import mp
  462. >>> A = mp.matrix([[1, -3 - 1j], [-3 + 1j, -2]])
  463. >>> E = mp.eighe(A, eigvals_only = True)
  464. >>> print(E)
  465. [-4.0]
  466. [ 3.0]
  467. >>> A = mp.matrix([[1, 2 + 5j], [2 - 5j, 3]])
  468. >>> E, Q = mp.eighe(A)
  469. >>> print(mp.chop(A * Q[:,0] - E[0] * Q[:,0]))
  470. [0.0]
  471. [0.0]
  472. see also: eigsy, eigh, eig
  473. """
  474. if not overwrite_a:
  475. A = A.copy()
  476. d = ctx.zeros(A.rows, 1)
  477. e = ctx.zeros(A.rows, 1)
  478. t = ctx.zeros(A.rows, 1)
  479. if eigvals_only:
  480. c_he_tridiag_0(ctx, A, d, e, t)
  481. tridiag_eigen(ctx, d, e, False)
  482. return d
  483. else:
  484. c_he_tridiag_0(ctx, A, d, e, t)
  485. B = ctx.eye(A.rows)
  486. tridiag_eigen(ctx, d, e, B)
  487. c_he_tridiag_2(ctx, A, t, B)
  488. return (d, B)
  489. @defun
  490. def eigh(ctx, A, eigvals_only = False, overwrite_a = False):
  491. """
  492. "eigh" is a unified interface for "eigsy" and "eighe". Depending on
  493. whether A is real or complex the appropriate function is called.
  494. This routine solves the (ordinary) eigenvalue problem for a real symmetric
  495. or complex hermitian square matrix A. Given A, an orthogonal (A real) or
  496. unitary (A complex) matrix Q is calculated which diagonalizes A:
  497. Q' A Q = diag(E) and Q Q' = Q' Q = 1
  498. Here diag(E) a is diagonal matrix whose diagonal is E.
  499. ' denotes the hermitian transpose (i.e. ordinary transposition and
  500. complex conjugation).
  501. The columns of Q are the eigenvectors of A and E contains the eigenvalues:
  502. A Q[:,i] = E[i] Q[:,i]
  503. input:
  504. A: a real or complex square matrix of format (n,n) which is symmetric
  505. (i.e. A[i,j]=A[j,i]) or hermitian (i.e. A[i,j]=conj(A[j,i])).
  506. eigvals_only: if true, calculates only the eigenvalues E.
  507. if false, calculates both eigenvectors and eigenvalues.
  508. overwrite_a: if true, allows modification of A which may improve
  509. performance. if false, A is not modified.
  510. output:
  511. E: vector of format (n). contains the eigenvalues of A in ascending order.
  512. Q: an orthogonal or unitary matrix of format (n,n). contains the
  513. eigenvectors of A as columns.
  514. return value:
  515. E if eigvals_only is true
  516. (E, Q) if eigvals_only is false
  517. example:
  518. >>> from mpmath import mp
  519. >>> A = mp.matrix([[3, 2], [2, 0]])
  520. >>> E = mp.eigh(A, eigvals_only = True)
  521. >>> print(E)
  522. [-1.0]
  523. [ 4.0]
  524. >>> A = mp.matrix([[1, 2], [2, 3]])
  525. >>> E, Q = mp.eigh(A)
  526. >>> print(mp.chop(A * Q[:,0] - E[0] * Q[:,0]))
  527. [0.0]
  528. [0.0]
  529. >>> A = mp.matrix([[1, 2 + 5j], [2 - 5j, 3]])
  530. >>> E, Q = mp.eigh(A)
  531. >>> print(mp.chop(A * Q[:,0] - E[0] * Q[:,0]))
  532. [0.0]
  533. [0.0]
  534. see also: eigsy, eighe, eig
  535. """
  536. iscomplex = any(type(x) is ctx.mpc for x in A)
  537. if iscomplex:
  538. return ctx.eighe(A, eigvals_only = eigvals_only, overwrite_a = overwrite_a)
  539. else:
  540. return ctx.eigsy(A, eigvals_only = eigvals_only, overwrite_a = overwrite_a)
  541. @defun
  542. def gauss_quadrature(ctx, n, qtype = "legendre", alpha = 0, beta = 0):
  543. """
  544. This routine calulates gaussian quadrature rules for different
  545. families of orthogonal polynomials. Let (a, b) be an interval,
  546. W(x) a positive weight function and n a positive integer.
  547. Then the purpose of this routine is to calculate pairs (x_k, w_k)
  548. for k=0, 1, 2, ... (n-1) which give
  549. int(W(x) * F(x), x = a..b) = sum(w_k * F(x_k),k = 0..(n-1))
  550. exact for all polynomials F(x) of degree (strictly) less than 2*n. For all
  551. integrable functions F(x) the sum is a (more or less) good approximation to
  552. the integral. The x_k are called nodes (which are the zeros of the
  553. related orthogonal polynomials) and the w_k are called the weights.
  554. parameters
  555. n (input) The degree of the quadrature rule, i.e. its number of
  556. nodes.
  557. qtype (input) The family of orthogonal polynmomials for which to
  558. compute the quadrature rule. See the list below.
  559. alpha (input) real number, used as parameter for some orthogonal
  560. polynomials
  561. beta (input) real number, used as parameter for some orthogonal
  562. polynomials.
  563. return value
  564. (X, W) a pair of two real arrays where x_k = X[k] and w_k = W[k].
  565. orthogonal polynomials:
  566. qtype polynomial
  567. ----- ----------
  568. "legendre" Legendre polynomials, W(x)=1 on the interval (-1, +1)
  569. "legendre01" shifted Legendre polynomials, W(x)=1 on the interval (0, +1)
  570. "hermite" Hermite polynomials, W(x)=exp(-x*x) on (-infinity,+infinity)
  571. "laguerre" Laguerre polynomials, W(x)=exp(-x) on (0,+infinity)
  572. "glaguerre" generalized Laguerre polynomials, W(x)=exp(-x)*x**alpha
  573. on (0, +infinity)
  574. "chebyshev1" Chebyshev polynomials of the first kind, W(x)=1/sqrt(1-x*x)
  575. on (-1, +1)
  576. "chebyshev2" Chebyshev polynomials of the second kind, W(x)=sqrt(1-x*x)
  577. on (-1, +1)
  578. "jacobi" Jacobi polynomials, W(x)=(1-x)**alpha * (1+x)**beta on (-1, +1)
  579. with alpha>-1 and beta>-1
  580. examples:
  581. >>> from mpmath import mp
  582. >>> f = lambda x: x**8 + 2 * x**6 - 3 * x**4 + 5 * x**2 - 7
  583. >>> X, W = mp.gauss_quadrature(5, "hermite")
  584. >>> A = mp.fdot([(f(x), w) for x, w in zip(X, W)])
  585. >>> B = mp.sqrt(mp.pi) * 57 / 16
  586. >>> C = mp.quad(lambda x: mp.exp(- x * x) * f(x), [-mp.inf, +mp.inf])
  587. >>> mp.nprint((mp.chop(A-B, tol = 1e-10), mp.chop(A-C, tol = 1e-10)))
  588. (0.0, 0.0)
  589. >>> f = lambda x: x**5 - 2 * x**4 + 3 * x**3 - 5 * x**2 + 7 * x - 11
  590. >>> X, W = mp.gauss_quadrature(3, "laguerre")
  591. >>> A = mp.fdot([(f(x), w) for x, w in zip(X, W)])
  592. >>> B = 76
  593. >>> C = mp.quad(lambda x: mp.exp(-x) * f(x), [0, +mp.inf])
  594. >>> mp.nprint(mp.chop(A-B, tol = 1e-10), mp.chop(A-C, tol = 1e-10))
  595. .0
  596. # orthogonality of the chebyshev polynomials:
  597. >>> f = lambda x: mp.chebyt(3, x) * mp.chebyt(2, x)
  598. >>> X, W = mp.gauss_quadrature(3, "chebyshev1")
  599. >>> A = mp.fdot([(f(x), w) for x, w in zip(X, W)])
  600. >>> print(mp.chop(A, tol = 1e-10))
  601. 0.0
  602. references:
  603. - golub and welsch, "calculations of gaussian quadrature rules", mathematics of
  604. computation 23, p. 221-230 (1969)
  605. - golub, "some modified matrix eigenvalue problems", siam review 15, p. 318-334 (1973)
  606. - stroud and secrest, "gaussian quadrature formulas", prentice-hall (1966)
  607. See also the routine gaussq.f in netlog.org or ACM Transactions on
  608. Mathematical Software algorithm 726.
  609. """
  610. d = ctx.zeros(n, 1)
  611. e = ctx.zeros(n, 1)
  612. z = ctx.zeros(1, n)
  613. z[0,0] = 1
  614. if qtype == "legendre":
  615. # legendre on the range -1 +1 , abramowitz, table 25.4, p.916
  616. w = 2
  617. for i in xrange(n):
  618. j = i + 1
  619. e[i] = ctx.sqrt(j * j / (4 * j * j - ctx.mpf(1)))
  620. elif qtype == "legendre01":
  621. # legendre shifted to 0 1 , abramowitz, table 25.8, p.921
  622. w = 1
  623. for i in xrange(n):
  624. d[i] = 1 / ctx.mpf(2)
  625. j = i + 1
  626. e[i] = ctx.sqrt(j * j / (16 * j * j - ctx.mpf(4)))
  627. elif qtype == "hermite":
  628. # hermite on the range -inf +inf , abramowitz, table 25.10,p.924
  629. w = ctx.sqrt(ctx.pi)
  630. for i in xrange(n):
  631. j = i + 1
  632. e[i] = ctx.sqrt(j / ctx.mpf(2))
  633. elif qtype == "laguerre":
  634. # laguerre on the range 0 +inf , abramowitz, table 25.9, p. 923
  635. w = 1
  636. for i in xrange(n):
  637. j = i + 1
  638. d[i] = 2 * j - 1
  639. e[i] = j
  640. elif qtype=="chebyshev1":
  641. # chebyshev polynimials of the first kind
  642. w = ctx.pi
  643. for i in xrange(n):
  644. e[i] = 1 / ctx.mpf(2)
  645. e[0] = ctx.sqrt(1 / ctx.mpf(2))
  646. elif qtype == "chebyshev2":
  647. # chebyshev polynimials of the second kind
  648. w = ctx.pi / 2
  649. for i in xrange(n):
  650. e[i] = 1 / ctx.mpf(2)
  651. elif qtype == "glaguerre":
  652. # generalized laguerre on the range 0 +inf
  653. w = ctx.gamma(1 + alpha)
  654. for i in xrange(n):
  655. j = i + 1
  656. d[i] = 2 * j - 1 + alpha
  657. e[i] = ctx.sqrt(j * (j + alpha))
  658. elif qtype == "jacobi":
  659. # jacobi polynomials
  660. alpha = ctx.mpf(alpha)
  661. beta = ctx.mpf(beta)
  662. ab = alpha + beta
  663. abi = ab + 2
  664. w = (2**(ab+1)) * ctx.gamma(alpha + 1) * ctx.gamma(beta + 1) / ctx.gamma(abi)
  665. d[0] = (beta - alpha) / abi
  666. e[0] = ctx.sqrt(4 * (1 + alpha) * (1 + beta) / ((abi + 1) * (abi * abi)))
  667. a2b2 = beta * beta - alpha * alpha
  668. for i in xrange(1, n):
  669. j = i + 1
  670. abi = 2 * j + ab
  671. d[i] = a2b2 / ((abi - 2) * abi)
  672. e[i] = ctx.sqrt(4 * j * (j + alpha) * (j + beta) * (j + ab) / ((abi * abi - 1) * abi * abi))
  673. elif isinstance(qtype, str):
  674. raise ValueError("unknown quadrature rule \"%s\"" % qtype)
  675. elif not isinstance(qtype, str):
  676. w = qtype(d, e)
  677. else:
  678. assert 0
  679. tridiag_eigen(ctx, d, e, z)
  680. for i in xrange(len(z)):
  681. z[i] *= z[i]
  682. z = z.transpose()
  683. return (d, w * z)
  684. ##################################################################################################
  685. ##################################################################################################
  686. ##################################################################################################
  687. def svd_r_raw(ctx, A, V = False, calc_u = False):
  688. """
  689. This routine computes the singular value decomposition of a matrix A.
  690. Given A, two orthogonal matrices U and V are calculated such that
  691. A = U S V
  692. where S is a suitable shaped matrix whose off-diagonal elements are zero.
  693. The diagonal elements of S are the singular values of A, i.e. the
  694. squareroots of the eigenvalues of A' A or A A'. Here ' denotes the transpose.
  695. Householder bidiagonalization and a variant of the QR algorithm is used.
  696. overview of the matrices :
  697. A : m*n A gets replaced by U
  698. U : m*n U replaces A. If n>m then only the first m*m block of U is
  699. non-zero. column-orthogonal: U' U = B
  700. here B is a n*n matrix whose first min(m,n) diagonal
  701. elements are 1 and all other elements are zero.
  702. S : n*n diagonal matrix, only the diagonal elements are stored in
  703. the array S. only the first min(m,n) diagonal elements are non-zero.
  704. V : n*n orthogonal: V V' = V' V = 1
  705. parameters:
  706. A (input/output) On input, A contains a real matrix of shape m*n.
  707. On output, if calc_u is true A contains the column-orthogonal
  708. matrix U; otherwise A is simply used as workspace and thus destroyed.
  709. V (input/output) if false, the matrix V is not calculated. otherwise
  710. V must be a matrix of shape n*n.
  711. calc_u (input) If true, the matrix U is calculated and replaces A.
  712. if false, U is not calculated and A is simply destroyed
  713. return value:
  714. S an array of length n containing the singular values of A sorted by
  715. decreasing magnitude. only the first min(m,n) elements are non-zero.
  716. This routine is a python translation of the fortran routine svd.f in the
  717. software library EISPACK (see netlib.org) which itself is based on the
  718. algol procedure svd described in:
  719. - num. math. 14, 403-420(1970) by golub and reinsch.
  720. - wilkinson/reinsch: handbook for auto. comp., vol ii-linear algebra, 134-151(1971).
  721. """
  722. m, n = A.rows, A.cols
  723. S = ctx.zeros(n, 1)
  724. # work is a temporary array of size n
  725. work = ctx.zeros(n, 1)
  726. g = scale = anorm = 0
  727. maxits = 3 * ctx.dps
  728. for i in xrange(n): # householder reduction to bidiagonal form
  729. work[i] = scale*g
  730. g = s = scale = 0
  731. if i < m:
  732. for k in xrange(i, m):
  733. scale += ctx.fabs(A[k,i])
  734. if scale != 0:
  735. for k in xrange(i, m):
  736. A[k,i] /= scale
  737. s += A[k,i] * A[k,i]
  738. f = A[i,i]
  739. g = -ctx.sqrt(s)
  740. if f < 0:
  741. g = -g
  742. h = f * g - s
  743. A[i,i] = f - g
  744. for j in xrange(i+1, n):
  745. s = 0
  746. for k in xrange(i, m):
  747. s += A[k,i] * A[k,j]
  748. f = s / h
  749. for k in xrange(i, m):
  750. A[k,j] += f * A[k,i]
  751. for k in xrange(i,m):
  752. A[k,i] *= scale
  753. S[i] = scale * g
  754. g = s = scale = 0
  755. if i < m and i != n - 1:
  756. for k in xrange(i+1, n):
  757. scale += ctx.fabs(A[i,k])
  758. if scale:
  759. for k in xrange(i+1, n):
  760. A[i,k] /= scale
  761. s += A[i,k] * A[i,k]
  762. f = A[i,i+1]
  763. g = -ctx.sqrt(s)
  764. if f < 0:
  765. g = -g
  766. h = f * g - s
  767. A[i,i+1] = f - g
  768. for k in xrange(i+1, n):
  769. work[k] = A[i,k] / h
  770. for j in xrange(i+1, m):
  771. s = 0
  772. for k in xrange(i+1, n):
  773. s += A[j,k] * A[i,k]
  774. for k in xrange(i+1, n):
  775. A[j,k] += s * work[k]
  776. for k in xrange(i+1, n):
  777. A[i,k] *= scale
  778. anorm = max(anorm, ctx.fabs(S[i]) + ctx.fabs(work[i]))
  779. if not isinstance(V, bool):
  780. for i in xrange(n-2, -1, -1): # accumulation of right hand transformations
  781. V[i+1,i+1] = 1
  782. if work[i+1] != 0:
  783. for j in xrange(i+1, n):
  784. V[i,j] = (A[i,j] / A[i,i+1]) / work[i+1]
  785. for j in xrange(i+1, n):
  786. s = 0
  787. for k in xrange(i+1, n):
  788. s += A[i,k] * V[j,k]
  789. for k in xrange(i+1, n):
  790. V[j,k] += s * V[i,k]
  791. for j in xrange(i+1, n):
  792. V[j,i] = V[i,j] = 0
  793. V[0,0] = 1
  794. if m<n : minnm = m
  795. else : minnm = n
  796. if calc_u:
  797. for i in xrange(minnm-1, -1, -1): # accumulation of left hand transformations
  798. g = S[i]
  799. for j in xrange(i+1, n):
  800. A[i,j] = 0
  801. if g != 0:
  802. g = 1 / g
  803. for j in xrange(i+1, n):
  804. s = 0
  805. for k in xrange(i+1, m):
  806. s += A[k,i] * A[k,j]
  807. f = (s / A[i,i]) * g
  808. for k in xrange(i, m):
  809. A[k,j] += f * A[k,i]
  810. for j in xrange(i, m):
  811. A[j,i] *= g
  812. else:
  813. for j in xrange(i, m):
  814. A[j,i] = 0
  815. A[i,i] += 1
  816. for k in xrange(n - 1, -1, -1):
  817. # diagonalization of the bidiagonal form:
  818. # loop over singular values, and over allowed itations
  819. its = 0
  820. while 1:
  821. its += 1
  822. flag = True
  823. for l in xrange(k, -1, -1):
  824. nm = l-1
  825. if ctx.fabs(work[l]) + anorm == anorm:
  826. flag = False
  827. break
  828. if ctx.fabs(S[nm]) + anorm == anorm:
  829. break
  830. if flag:
  831. c = 0
  832. s = 1
  833. for i in xrange(l, k + 1):
  834. f = s * work[i]
  835. work[i] *= c
  836. if ctx.fabs(f) + anorm == anorm:
  837. break
  838. g = S[i]
  839. h = ctx.hypot(f, g)
  840. S[i] = h
  841. h = 1 / h
  842. c = g * h
  843. s = - f * h
  844. if calc_u:
  845. for j in xrange(m):
  846. y = A[j,nm]
  847. z = A[j,i]
  848. A[j,nm] = y * c + z * s
  849. A[j,i] = z * c - y * s
  850. z = S[k]
  851. if l == k: # convergence
  852. if z < 0: # singular value is made nonnegative
  853. S[k] = -z
  854. if not isinstance(V, bool):
  855. for j in xrange(n):
  856. V[k,j] = -V[k,j]
  857. break
  858. if its >= maxits:
  859. raise RuntimeError("svd: no convergence to an eigenvalue after %d iterations" % its)
  860. x = S[l] # shift from bottom 2 by 2 minor
  861. nm = k-1
  862. y = S[nm]
  863. g = work[nm]
  864. h = work[k]
  865. f = ((y - z) * (y + z) + (g - h) * (g + h))/(2 * h * y)
  866. g = ctx.hypot(f, 1)
  867. if f >= 0: f = ((x - z) * (x + z) + h * ((y / (f + g)) - h)) / x
  868. else: f = ((x - z) * (x + z) + h * ((y / (f - g)) - h)) / x
  869. c = s = 1 # next qt transformation
  870. for j in xrange(l, nm + 1):
  871. g = work[j+1]
  872. y = S[j+1]
  873. h = s * g
  874. g = c * g
  875. z = ctx.hypot(f, h)
  876. work[j] = z
  877. c = f / z
  878. s = h / z
  879. f = x * c + g * s
  880. g = g * c - x * s
  881. h = y * s
  882. y *= c
  883. if not isinstance(V, bool):
  884. for jj in xrange(n):
  885. x = V[j ,jj]
  886. z = V[j+1,jj]
  887. V[j ,jj]= x * c + z * s
  888. V[j+1 ,jj]= z * c - x * s
  889. z = ctx.hypot(f, h)
  890. S[j] = z
  891. if z != 0: # rotation can be arbitray if z=0
  892. z = 1 / z
  893. c = f * z
  894. s = h * z
  895. f = c * g + s * y
  896. x = c * y - s * g
  897. if calc_u:
  898. for jj in xrange(m):
  899. y = A[jj,j ]
  900. z = A[jj,j+1]
  901. A[jj,j ] = y * c + z * s
  902. A[jj,j+1 ] = z * c - y * s
  903. work[l] = 0
  904. work[k] = f
  905. S[k] = x
  906. ##########################
  907. # Sort singular values into decreasing order (bubble-sort)
  908. for i in xrange(n):
  909. imax = i
  910. s = ctx.fabs(S[i]) # s is the current maximal element
  911. for j in xrange(i + 1, n):
  912. c = ctx.fabs(S[j])
  913. if c > s:
  914. s = c
  915. imax = j
  916. if imax != i:
  917. # swap singular values
  918. z = S[i]
  919. S[i] = S[imax]
  920. S[imax] = z
  921. if calc_u:
  922. for j in xrange(m):
  923. z = A[j,i]
  924. A[j,i] = A[j,imax]
  925. A[j,imax] = z
  926. if not isinstance(V, bool):
  927. for j in xrange(n):
  928. z = V[i,j]
  929. V[i,j] = V[imax,j]
  930. V[imax,j] = z
  931. return S
  932. #######################
  933. def svd_c_raw(ctx, A, V = False, calc_u = False):
  934. """
  935. This routine computes the singular value decomposition of a matrix A.
  936. Given A, two unitary matrices U and V are calculated such that
  937. A = U S V
  938. where S is a suitable shaped matrix whose off-diagonal elements are zero.
  939. The diagonal elements of S are the singular values of A, i.e. the
  940. squareroots of the eigenvalues of A' A or A A'. Here ' denotes the hermitian
  941. transpose (i.e. transposition and conjugation). Householder bidiagonalization
  942. and a variant of the QR algorithm is used.
  943. overview of the matrices :
  944. A : m*n A gets replaced by U
  945. U : m*n U replaces A. If n>m then only the first m*m block of U is
  946. non-zero. column-unitary: U' U = B
  947. here B is a n*n matrix whose first min(m,n) diagonal
  948. elements are 1 and all other elements are zero.
  949. S : n*n diagonal matrix, only the diagonal elements are stored in
  950. the array S. only the first min(m,n) diagonal elements are non-zero.
  951. V : n*n unitary: V V' = V' V = 1
  952. parameters:
  953. A (input/output) On input, A contains a complex matrix of shape m*n.
  954. On output, if calc_u is true A contains the column-unitary
  955. matrix U; otherwise A is simply used as workspace and thus destroyed.
  956. V (input/output) if false, the matrix V is not calculated. otherwise
  957. V must be a matrix of shape n*n.
  958. calc_u (input) If true, the matrix U is calculated and replaces A.
  959. if false, U is not calculated and A is simply destroyed
  960. return value:
  961. S an array of length n containing the singular values of A sorted by
  962. decreasing magnitude. only the first min(m,n) elements are non-zero.
  963. This routine is a python translation of the fortran routine svd.f in the
  964. software library EISPACK (see netlib.org) which itself is based on the
  965. algol procedure svd described in:
  966. - num. math. 14, 403-420(1970) by golub and reinsch.
  967. - wilkinson/reinsch: handbook for auto. comp., vol ii-linear algebra, 134-151(1971).
  968. """
  969. m, n = A.rows, A.cols
  970. S = ctx.zeros(n, 1)
  971. # work is a temporary array of size n
  972. work = ctx.zeros(n, 1)
  973. lbeta = ctx.zeros(n, 1)
  974. rbeta = ctx.zeros(n, 1)
  975. dwork = ctx.zeros(n, 1)
  976. g = scale = anorm = 0
  977. maxits = 3 * ctx.dps
  978. for i in xrange(n): # householder reduction to bidiagonal form
  979. dwork[i] = scale * g # dwork are the side-diagonal elements
  980. g = s = scale = 0
  981. if i < m:
  982. for k in xrange(i, m):
  983. scale += ctx.fabs(ctx.re(A[k,i])) + ctx.fabs(ctx.im(A[k,i]))
  984. if scale != 0:
  985. for k in xrange(i, m):
  986. A[k,i] /= scale
  987. ar = ctx.re(A[k,i])
  988. ai = ctx.im(A[k,i])
  989. s += ar * ar + ai * ai
  990. f = A[i,i]
  991. g = -ctx.sqrt(s)
  992. if ctx.re(f) < 0:
  993. beta = -g - ctx.conj(f)
  994. g = -g
  995. else:
  996. beta = -g + ctx.conj(f)
  997. beta /= ctx.conj(beta)
  998. beta += 1
  999. h = 2 * (ctx.re(f) * g - s)
  1000. A[i,i] = f - g
  1001. beta /= h
  1002. lbeta[i] = (beta / scale) / scale
  1003. for j in xrange(i+1, n):
  1004. s = 0
  1005. for k in xrange(i, m):
  1006. s += ctx.conj(A[k,i]) * A[k,j]
  1007. f = beta * s
  1008. for k in xrange(i, m):
  1009. A[k,j] += f * A[k,i]
  1010. for k in xrange(i, m):
  1011. A[k,i] *= scale
  1012. S[i] = scale * g # S are the diagonal elements
  1013. g = s = scale = 0
  1014. if i < m and i != n - 1:
  1015. for k in xrange(i+1, n):
  1016. scale += ctx.fabs(ctx.re(A[i,k])) + ctx.fabs(ctx.im(A[i,k]))
  1017. if scale:
  1018. for k in xrange(i+1, n):
  1019. A[i,k] /= scale
  1020. ar = ctx.re(A[i,k])
  1021. ai = ctx.im(A[i,k])
  1022. s += ar * ar + ai * ai
  1023. f = A[i,i+1]
  1024. g = -ctx.sqrt(s)
  1025. if ctx.re(f) < 0:
  1026. beta = -g - ctx.conj(f)
  1027. g = -g
  1028. else:
  1029. beta = -g + ctx.conj(f)
  1030. beta /= ctx.conj(beta)
  1031. beta += 1
  1032. h = 2 * (ctx.re(f) * g - s)
  1033. A[i,i+1] = f - g
  1034. beta /= h
  1035. rbeta[i] = (beta / scale) / scale
  1036. for k in xrange(i+1, n):
  1037. work[k] = A[i, k]
  1038. for j in xrange(i+1, m):
  1039. s = 0
  1040. for k in xrange(i+1, n):
  1041. s += ctx.conj(A[i,k]) * A[j,k]
  1042. f = s * beta
  1043. for k in xrange(i+1,n):
  1044. A[j,k] += f * work[k]
  1045. for k in xrange(i+1, n):
  1046. A[i,k] *= scale
  1047. anorm = max(anorm,ctx.fabs(S[i]) + ctx.fabs(dwork[i]))
  1048. if not isinstance(V, bool):
  1049. for i in xrange(n-2, -1, -1): # accumulation of right hand transformations
  1050. V[i+1,i+1] = 1
  1051. if dwork[i+1] != 0:
  1052. f = ctx.conj(rbeta[i])
  1053. for j in xrange(i+1, n):
  1054. V[i,j] = A[i,j] * f
  1055. for j in xrange(i+1, n):
  1056. s = 0
  1057. for k in xrange(i+1, n):
  1058. s += ctx.conj(A[i,k]) * V[j,k]
  1059. for k in xrange(i+1, n):
  1060. V[j,k] += s * V[i,k]
  1061. for j in xrange(i+1,n):
  1062. V[j,i] = V[i,j] = 0
  1063. V[0,0] = 1
  1064. if m < n : minnm = m
  1065. else : minnm = n
  1066. if calc_u:
  1067. for i in xrange(minnm-1, -1, -1): # accumulation of left hand transformations
  1068. g = S[i]
  1069. for j in xrange(i+1, n):
  1070. A[i,j] = 0
  1071. if g != 0:
  1072. g = 1 / g
  1073. for j in xrange(i+1, n):
  1074. s = 0
  1075. for k in xrange(i+1, m):
  1076. s += ctx.conj(A[k,i]) * A[k,j]
  1077. f = s * ctx.conj(lbeta[i])
  1078. for k in xrange(i, m):
  1079. A[k,j] += f * A[k,i]
  1080. for j in xrange(i, m):
  1081. A[j,i] *= g
  1082. else:
  1083. for j in xrange(i, m):
  1084. A[j,i] = 0
  1085. A[i,i] += 1
  1086. for k in xrange(n-1, -1, -1):
  1087. # diagonalization of the bidiagonal form:
  1088. # loop over singular values, and over allowed itations
  1089. its = 0
  1090. while 1:
  1091. its += 1
  1092. flag = True
  1093. for l in xrange(k, -1, -1):
  1094. nm = l - 1
  1095. if ctx.fabs(dwork[l]) + anorm == anorm:
  1096. flag = False
  1097. break
  1098. if ctx.fabs(S[nm]) + anorm == anorm:
  1099. break
  1100. if flag:
  1101. c = 0
  1102. s = 1
  1103. for i in xrange(l, k+1):
  1104. f = s * dwork[i]
  1105. dwork[i] *= c
  1106. if ctx.fabs(f) + anorm == anorm:
  1107. break
  1108. g = S[i]
  1109. h = ctx.hypot(f, g)
  1110. S[i] = h
  1111. h = 1 / h
  1112. c = g * h
  1113. s = -f * h
  1114. if calc_u:
  1115. for j in xrange(m):
  1116. y = A[j,nm]
  1117. z = A[j,i]
  1118. A[j,nm]= y * c + z * s
  1119. A[j,i] = z * c - y * s
  1120. z = S[k]
  1121. if l == k: # convergence
  1122. if z < 0: # singular value is made nonnegative
  1123. S[k] = -z
  1124. if not isinstance(V, bool):
  1125. for j in xrange(n):
  1126. V[k,j] = -V[k,j]
  1127. break
  1128. if its >= maxits:
  1129. raise RuntimeError("svd: no convergence to an eigenvalue after %d iterations" % its)
  1130. x = S[l] # shift from bottom 2 by 2 minor
  1131. nm = k-1
  1132. y = S[nm]
  1133. g = dwork[nm]
  1134. h = dwork[k]
  1135. f = ((y - z) * (y + z) + (g - h) * (g + h)) / (2 * h * y)
  1136. g = ctx.hypot(f, 1)
  1137. if f >=0: f = (( x - z) *( x + z) + h *((y / (f + g)) - h)) / x
  1138. else: f = (( x - z) *( x + z) + h *((y / (f - g)) - h)) / x
  1139. c = s = 1 # next qt transformation
  1140. for j in xrange(l, nm + 1):
  1141. g = dwork[j+1]
  1142. y = S[j+1]
  1143. h = s * g
  1144. g = c * g
  1145. z = ctx.hypot(f, h)
  1146. dwork[j] = z
  1147. c = f / z
  1148. s = h / z
  1149. f = x * c + g * s
  1150. g = g * c - x * s
  1151. h = y * s
  1152. y *= c
  1153. if not isinstance(V, bool):
  1154. for jj in xrange(n):
  1155. x = V[j ,jj]
  1156. z = V[j+1,jj]
  1157. V[j ,jj]= x * c + z * s
  1158. V[j+1,jj ]= z * c - x * s
  1159. z = ctx.hypot(f, h)
  1160. S[j] = z
  1161. if z != 0: # rotation can be arbitray if z=0
  1162. z = 1 / z
  1163. c = f * z
  1164. s = h * z
  1165. f = c * g + s * y
  1166. x = c * y - s * g
  1167. if calc_u:
  1168. for jj in xrange(m):
  1169. y = A[jj,j ]
  1170. z = A[jj,j+1]
  1171. A[jj,j ]= y * c + z * s
  1172. A[jj,j+1 ]= z * c - y * s
  1173. dwork[l] = 0
  1174. dwork[k] = f
  1175. S[k] = x
  1176. ##########################
  1177. # Sort singular values into decreasing order (bubble-sort)
  1178. for i in xrange(n):
  1179. imax = i
  1180. s = ctx.fabs(S[i]) # s is the current maximal element
  1181. for j in xrange(i + 1, n):
  1182. c = ctx.fabs(S[j])
  1183. if c > s:
  1184. s = c
  1185. imax = j
  1186. if imax != i:
  1187. # swap singular values
  1188. z = S[i]
  1189. S[i] = S[imax]
  1190. S[imax] = z
  1191. if calc_u:
  1192. for j in xrange(m):
  1193. z = A[j,i]
  1194. A[j,i] = A[j,imax]
  1195. A[j,imax] = z
  1196. if not isinstance(V, bool):
  1197. for j in xrange(n):
  1198. z = V[i,j]
  1199. V[i,j] = V[imax,j]
  1200. V[imax,j] = z
  1201. return S
  1202. ##################################################################################################
  1203. @defun
  1204. def svd_r(ctx, A, full_matrices = False, compute_uv = True, overwrite_a = False):
  1205. """
  1206. This routine computes the singular value decomposition of a matrix A.
  1207. Given A, two orthogonal matrices U and V are calculated such that
  1208. A = U S V and U' U = 1 and V V' = 1
  1209. where S is a suitable shaped matrix whose off-diagonal elements are zero.
  1210. Here ' denotes the transpose. The diagonal elements of S are the singular
  1211. values of A, i.e. the squareroots of the eigenvalues of A' A or A A'.
  1212. input:
  1213. A : a real matrix of shape (m, n)
  1214. full_matrices : if true, U and V are of shape (m, m) and (n, n).
  1215. if false, U and V are of shape (m, min(m, n)) and (min(m, n), n).
  1216. compute_uv : if true, U and V are calculated. if false, only S is calculated.
  1217. overwrite_a : if true, allows modification of A which may improve
  1218. performance. if false, A is not modified.
  1219. output:
  1220. U : an orthogonal matrix: U' U = 1. if full_matrices is true, U is of
  1221. shape (m, m). ortherwise it is of shape (m, min(m, n)).
  1222. S : an array of length min(m, n) containing the singular values of A sorted by
  1223. decreasing magnitude.
  1224. V : an orthogonal matrix: V V' = 1. if full_matrices is true, V is of
  1225. shape (n, n). ortherwise it is of shape (min(m, n), n).
  1226. return value:
  1227. S if compute_uv is false
  1228. (U, S, V) if compute_uv is true
  1229. overview of the matrices:
  1230. full_matrices true:
  1231. A : m*n
  1232. U : m*m U' U = 1
  1233. S as matrix : m*n
  1234. V : n*n V V' = 1
  1235. full_matrices false:
  1236. A : m*n
  1237. U : m*min(n,m) U' U = 1
  1238. S as matrix : min(m,n)*min(m,n)
  1239. V : min(m,n)*n V V' = 1
  1240. examples:
  1241. >>> from mpmath import mp
  1242. >>> A = mp.matrix([[2, -2, -1], [3, 4, -2], [-2, -2, 0]])
  1243. >>> S = mp.svd_r(A, compute_uv = False)
  1244. >>> print(S)
  1245. [6.0]
  1246. [3.0]
  1247. [1.0]
  1248. >>> U, S, V = mp.svd_r(A)
  1249. >>> print(mp.chop(A - U * mp.diag(S) * V))
  1250. [0.0 0.0 0.0]
  1251. [0.0 0.0 0.0]
  1252. [0.0 0.0 0.0]
  1253. see also: svd, svd_c
  1254. """
  1255. m, n = A.rows, A.cols
  1256. if not compute_uv:
  1257. if not overwrite_a:
  1258. A = A.copy()
  1259. S = svd_r_raw(ctx, A, V = False, calc_u = False)
  1260. S = S[:min(m,n)]
  1261. return S
  1262. if full_matrices and n < m:
  1263. V = ctx.zeros(m, m)
  1264. A0 = ctx.zeros(m, m)
  1265. A0[:,:n] = A
  1266. S = svd_r_raw(ctx, A0, V, calc_u = True)
  1267. S = S[:n]
  1268. V = V[:n,:n]
  1269. return (A0, S, V)
  1270. else:
  1271. if not overwrite_a:
  1272. A = A.copy()
  1273. V = ctx.zeros(n, n)
  1274. S = svd_r_raw(ctx, A, V, calc_u = True)
  1275. if n > m:
  1276. if full_matrices == False:
  1277. V = V[:m,:]
  1278. S = S[:m]
  1279. A = A[:,:m]
  1280. return (A, S, V)
  1281. ##############################
  1282. @defun
  1283. def svd_c(ctx, A, full_matrices = False, compute_uv = True, overwrite_a = False):
  1284. """
  1285. This routine computes the singular value decomposition of a matrix A.
  1286. Given A, two unitary matrices U and V are calculated such that
  1287. A = U S V and U' U = 1 and V V' = 1
  1288. where S is a suitable shaped matrix whose off-diagonal elements are zero.
  1289. Here ' denotes the hermitian transpose (i.e. transposition and complex
  1290. conjugation). The diagonal elements of S are the singular values of A,
  1291. i.e. the squareroots of the eigenvalues of A' A or A A'.
  1292. input:
  1293. A : a complex matrix of shape (m, n)
  1294. full_matrices : if true, U and V are of shape (m, m) and (n, n).
  1295. if false, U and V are of shape (m, min(m, n)) and (min(m, n), n).
  1296. compute_uv : if true, U and V are calculated. if false, only S is calculated.
  1297. overwrite_a : if true, allows modification of A which may improve
  1298. performance. if false, A is not modified.
  1299. output:
  1300. U : an unitary matrix: U' U = 1. if full_matrices is true, U is of
  1301. shape (m, m). ortherwise it is of shape (m, min(m, n)).
  1302. S : an array of length min(m, n) containing the singular values of A sorted by
  1303. decreasing magnitude.
  1304. V : an unitary matrix: V V' = 1. if full_matrices is true, V is of
  1305. shape (n, n). ortherwise it is of shape (min(m, n), n).
  1306. return value:
  1307. S if compute_uv is false
  1308. (U, S, V) if compute_uv is true
  1309. overview of the matrices:
  1310. full_matrices true:
  1311. A : m*n
  1312. U : m*m U' U = 1
  1313. S as matrix : m*n
  1314. V : n*n V V' = 1
  1315. full_matrices false:
  1316. A : m*n
  1317. U : m*min(n,m) U' U = 1
  1318. S as matrix : min(m,n)*min(m,n)
  1319. V : min(m,n)*n V V' = 1
  1320. example:
  1321. >>> from mpmath import mp
  1322. >>> A = mp.matrix([[-2j, -1-3j, -2+2j], [2-2j, -1-3j, 1], [-3+1j,-2j,0]])
  1323. >>> S = mp.svd_c(A, compute_uv = False)
  1324. >>> print(mp.chop(S - mp.matrix([mp.sqrt(34), mp.sqrt(15), mp.sqrt(6)])))
  1325. [0.0]
  1326. [0.0]
  1327. [0.0]
  1328. >>> U, S, V = mp.svd_c(A)
  1329. >>> print(mp.chop(A - U * mp.diag(S) * V))
  1330. [0.0 0.0 0.0]
  1331. [0.0 0.0 0.0]
  1332. [0.0 0.0 0.0]
  1333. see also: svd, svd_r
  1334. """
  1335. m, n = A.rows, A.cols
  1336. if not compute_uv:
  1337. if not overwrite_a:
  1338. A = A.copy()
  1339. S = svd_c_raw(ctx, A, V = False, calc_u = False)
  1340. S = S[:min(m,n)]
  1341. return S
  1342. if full_matrices and n < m:
  1343. V = ctx.zeros(m, m)
  1344. A0 = ctx.zeros(m, m)
  1345. A0[:,:n] = A
  1346. S = svd_c_raw(ctx, A0, V, calc_u = True)
  1347. S = S[:n]
  1348. V = V[:n,:n]
  1349. return (A0, S, V)
  1350. else:
  1351. if not overwrite_a:
  1352. A = A.copy()
  1353. V = ctx.zeros(n, n)
  1354. S = svd_c_raw(ctx, A, V, calc_u = True)
  1355. if n > m:
  1356. if full_matrices == False:
  1357. V = V[:m,:]
  1358. S = S[:m]
  1359. A = A[:,:m]
  1360. return (A, S, V)
  1361. @defun
  1362. def svd(ctx, A, full_matrices = False, compute_uv = True, overwrite_a = False):
  1363. """
  1364. "svd" is a unified interface for "svd_r" and "svd_c". Depending on
  1365. whether A is real or complex the appropriate function is called.
  1366. This routine computes the singular value decomposition of a matrix A.
  1367. Given A, two orthogonal (A real) or unitary (A complex) matrices U and V
  1368. are calculated such that
  1369. A = U S V and U' U = 1 and V V' = 1
  1370. where S is a suitable shaped matrix whose off-diagonal elements are zero.
  1371. Here ' denotes the hermitian transpose (i.e. transposition and complex
  1372. conjugation). The diagonal elements of S are the singular values of A,
  1373. i.e. the squareroots of the eigenvalues of A' A or A A'.
  1374. input:
  1375. A : a real or complex matrix of shape (m, n)
  1376. full_matrices : if true, U and V are of shape (m, m) and (n, n).
  1377. if false, U and V are of shape (m, min(m, n)) and (min(m, n), n).
  1378. compute_uv : if true, U and V are calculated. if false, only S is calculated.
  1379. overwrite_a : if true, allows modification of A which may improve
  1380. performance. if false, A is not modified.
  1381. output:
  1382. U : an orthogonal or unitary matrix: U' U = 1. if full_matrices is true, U is of
  1383. shape (m, m). ortherwise it is of shape (m, min(m, n)).
  1384. S : an array of length min(m, n) containing the singular values of A sorted by
  1385. decreasing magnitude.
  1386. V : an orthogonal or unitary matrix: V V' = 1. if full_matrices is true, V is of
  1387. shape (n, n). ortherwise it is of shape (min(m, n), n).
  1388. return value:
  1389. S if compute_uv is false
  1390. (U, S, V) if compute_uv is true
  1391. overview of the matrices:
  1392. full_matrices true:
  1393. A : m*n
  1394. U : m*m U' U = 1
  1395. S as matrix : m*n
  1396. V : n*n V V' = 1
  1397. full_matrices false:
  1398. A : m*n
  1399. U : m*min(n,m) U' U = 1
  1400. S as matrix : min(m,n)*min(m,n)
  1401. V : min(m,n)*n V V' = 1
  1402. examples:
  1403. >>> from mpmath import mp
  1404. >>> A = mp.matrix([[2, -2, -1], [3, 4, -2], [-2, -2, 0]])
  1405. >>> S = mp.svd(A, compute_uv = False)
  1406. >>> print(S)
  1407. [6.0]
  1408. [3.0]
  1409. [1.0]
  1410. >>> U, S, V = mp.svd(A)
  1411. >>> print(mp.chop(A - U * mp.diag(S) * V))
  1412. [0.0 0.0 0.0]
  1413. [0.0 0.0 0.0]
  1414. [0.0 0.0 0.0]
  1415. see also: svd_r, svd_c
  1416. """
  1417. iscomplex = any(type(x) is ctx.mpc for x in A)
  1418. if iscomplex:
  1419. return ctx.svd_c(A, full_matrices = full_matrices, compute_uv = compute_uv, overwrite_a = overwrite_a)
  1420. else:
  1421. return ctx.svd_r(A, full_matrices = full_matrices, compute_uv = compute_uv, overwrite_a = overwrite_a)