arpack.py 72 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892
  1. """
  2. Find a few eigenvectors and eigenvalues of a matrix.
  3. Uses ARPACK: http://www.caam.rice.edu/software/ARPACK/
  4. """
  5. # Wrapper implementation notes
  6. #
  7. # ARPACK Entry Points
  8. # -------------------
  9. # The entry points to ARPACK are
  10. # - (s,d)seupd : single and double precision symmetric matrix
  11. # - (s,d,c,z)neupd: single,double,complex,double complex general matrix
  12. # This wrapper puts the *neupd (general matrix) interfaces in eigs()
  13. # and the *seupd (symmetric matrix) in eigsh().
  14. # There is no Hermetian complex/double complex interface.
  15. # To find eigenvalues of a Hermetian matrix you
  16. # must use eigs() and not eigsh()
  17. # It might be desirable to handle the Hermetian case differently
  18. # and, for example, return real eigenvalues.
  19. # Number of eigenvalues returned and complex eigenvalues
  20. # ------------------------------------------------------
  21. # The ARPACK nonsymmetric real and double interface (s,d)naupd return
  22. # eigenvalues and eigenvectors in real (float,double) arrays.
  23. # Since the eigenvalues and eigenvectors are, in general, complex
  24. # ARPACK puts the real and imaginary parts in consecutive entries
  25. # in real-valued arrays. This wrapper puts the real entries
  26. # into complex data types and attempts to return the requested eigenvalues
  27. # and eigenvectors.
  28. # Solver modes
  29. # ------------
  30. # ARPACK and handle shifted and shift-inverse computations
  31. # for eigenvalues by providing a shift (sigma) and a solver.
  32. from __future__ import division, print_function, absolute_import
  33. __docformat__ = "restructuredtext en"
  34. __all__ = ['eigs', 'eigsh', 'svds', 'ArpackError', 'ArpackNoConvergence']
  35. from . import _arpack
  36. import numpy as np
  37. import warnings
  38. from scipy.sparse.linalg.interface import aslinearoperator, LinearOperator
  39. from scipy.sparse import eye, issparse, isspmatrix, isspmatrix_csr
  40. from scipy.linalg import eig, eigh, lu_factor, lu_solve
  41. from scipy.sparse.sputils import isdense
  42. from scipy.sparse.linalg import gmres, splu
  43. from scipy._lib._util import _aligned_zeros
  44. from scipy._lib._threadsafety import ReentrancyLock
  45. _type_conv = {'f': 's', 'd': 'd', 'F': 'c', 'D': 'z'}
  46. _ndigits = {'f': 5, 'd': 12, 'F': 5, 'D': 12}
  47. DNAUPD_ERRORS = {
  48. 0: "Normal exit.",
  49. 1: "Maximum number of iterations taken. "
  50. "All possible eigenvalues of OP has been found. IPARAM(5) "
  51. "returns the number of wanted converged Ritz values.",
  52. 2: "No longer an informational error. Deprecated starting "
  53. "with release 2 of ARPACK.",
  54. 3: "No shifts could be applied during a cycle of the "
  55. "Implicitly restarted Arnoldi iteration. One possibility "
  56. "is to increase the size of NCV relative to NEV. ",
  57. -1: "N must be positive.",
  58. -2: "NEV must be positive.",
  59. -3: "NCV-NEV >= 2 and less than or equal to N.",
  60. -4: "The maximum number of Arnoldi update iterations allowed "
  61. "must be greater than zero.",
  62. -5: " WHICH must be one of 'LM', 'SM', 'LR', 'SR', 'LI', 'SI'",
  63. -6: "BMAT must be one of 'I' or 'G'.",
  64. -7: "Length of private work array WORKL is not sufficient.",
  65. -8: "Error return from LAPACK eigenvalue calculation;",
  66. -9: "Starting vector is zero.",
  67. -10: "IPARAM(7) must be 1,2,3,4.",
  68. -11: "IPARAM(7) = 1 and BMAT = 'G' are incompatible.",
  69. -12: "IPARAM(1) must be equal to 0 or 1.",
  70. -13: "NEV and WHICH = 'BE' are incompatible.",
  71. -9999: "Could not build an Arnoldi factorization. "
  72. "IPARAM(5) returns the size of the current Arnoldi "
  73. "factorization. The user is advised to check that "
  74. "enough workspace and array storage has been allocated."
  75. }
  76. SNAUPD_ERRORS = DNAUPD_ERRORS
  77. ZNAUPD_ERRORS = DNAUPD_ERRORS.copy()
  78. ZNAUPD_ERRORS[-10] = "IPARAM(7) must be 1,2,3."
  79. CNAUPD_ERRORS = ZNAUPD_ERRORS
  80. DSAUPD_ERRORS = {
  81. 0: "Normal exit.",
  82. 1: "Maximum number of iterations taken. "
  83. "All possible eigenvalues of OP has been found.",
  84. 2: "No longer an informational error. Deprecated starting with "
  85. "release 2 of ARPACK.",
  86. 3: "No shifts could be applied during a cycle of the Implicitly "
  87. "restarted Arnoldi iteration. One possibility is to increase "
  88. "the size of NCV relative to NEV. ",
  89. -1: "N must be positive.",
  90. -2: "NEV must be positive.",
  91. -3: "NCV must be greater than NEV and less than or equal to N.",
  92. -4: "The maximum number of Arnoldi update iterations allowed "
  93. "must be greater than zero.",
  94. -5: "WHICH must be one of 'LM', 'SM', 'LA', 'SA' or 'BE'.",
  95. -6: "BMAT must be one of 'I' or 'G'.",
  96. -7: "Length of private work array WORKL is not sufficient.",
  97. -8: "Error return from trid. eigenvalue calculation; "
  98. "Informational error from LAPACK routine dsteqr .",
  99. -9: "Starting vector is zero.",
  100. -10: "IPARAM(7) must be 1,2,3,4,5.",
  101. -11: "IPARAM(7) = 1 and BMAT = 'G' are incompatible.",
  102. -12: "IPARAM(1) must be equal to 0 or 1.",
  103. -13: "NEV and WHICH = 'BE' are incompatible. ",
  104. -9999: "Could not build an Arnoldi factorization. "
  105. "IPARAM(5) returns the size of the current Arnoldi "
  106. "factorization. The user is advised to check that "
  107. "enough workspace and array storage has been allocated.",
  108. }
  109. SSAUPD_ERRORS = DSAUPD_ERRORS
  110. DNEUPD_ERRORS = {
  111. 0: "Normal exit.",
  112. 1: "The Schur form computed by LAPACK routine dlahqr "
  113. "could not be reordered by LAPACK routine dtrsen. "
  114. "Re-enter subroutine dneupd with IPARAM(5)NCV and "
  115. "increase the size of the arrays DR and DI to have "
  116. "dimension at least dimension NCV and allocate at least NCV "
  117. "columns for Z. NOTE: Not necessary if Z and V share "
  118. "the same space. Please notify the authors if this error"
  119. "occurs.",
  120. -1: "N must be positive.",
  121. -2: "NEV must be positive.",
  122. -3: "NCV-NEV >= 2 and less than or equal to N.",
  123. -5: "WHICH must be one of 'LM', 'SM', 'LR', 'SR', 'LI', 'SI'",
  124. -6: "BMAT must be one of 'I' or 'G'.",
  125. -7: "Length of private work WORKL array is not sufficient.",
  126. -8: "Error return from calculation of a real Schur form. "
  127. "Informational error from LAPACK routine dlahqr .",
  128. -9: "Error return from calculation of eigenvectors. "
  129. "Informational error from LAPACK routine dtrevc.",
  130. -10: "IPARAM(7) must be 1,2,3,4.",
  131. -11: "IPARAM(7) = 1 and BMAT = 'G' are incompatible.",
  132. -12: "HOWMNY = 'S' not yet implemented",
  133. -13: "HOWMNY must be one of 'A' or 'P' if RVEC = .true.",
  134. -14: "DNAUPD did not find any eigenvalues to sufficient "
  135. "accuracy.",
  136. -15: "DNEUPD got a different count of the number of converged "
  137. "Ritz values than DNAUPD got. This indicates the user "
  138. "probably made an error in passing data from DNAUPD to "
  139. "DNEUPD or that the data was modified before entering "
  140. "DNEUPD",
  141. }
  142. SNEUPD_ERRORS = DNEUPD_ERRORS.copy()
  143. SNEUPD_ERRORS[1] = ("The Schur form computed by LAPACK routine slahqr "
  144. "could not be reordered by LAPACK routine strsen . "
  145. "Re-enter subroutine dneupd with IPARAM(5)=NCV and "
  146. "increase the size of the arrays DR and DI to have "
  147. "dimension at least dimension NCV and allocate at least "
  148. "NCV columns for Z. NOTE: Not necessary if Z and V share "
  149. "the same space. Please notify the authors if this error "
  150. "occurs.")
  151. SNEUPD_ERRORS[-14] = ("SNAUPD did not find any eigenvalues to sufficient "
  152. "accuracy.")
  153. SNEUPD_ERRORS[-15] = ("SNEUPD got a different count of the number of "
  154. "converged Ritz values than SNAUPD got. This indicates "
  155. "the user probably made an error in passing data from "
  156. "SNAUPD to SNEUPD or that the data was modified before "
  157. "entering SNEUPD")
  158. ZNEUPD_ERRORS = {0: "Normal exit.",
  159. 1: "The Schur form computed by LAPACK routine csheqr "
  160. "could not be reordered by LAPACK routine ztrsen. "
  161. "Re-enter subroutine zneupd with IPARAM(5)=NCV and "
  162. "increase the size of the array D to have "
  163. "dimension at least dimension NCV and allocate at least "
  164. "NCV columns for Z. NOTE: Not necessary if Z and V share "
  165. "the same space. Please notify the authors if this error "
  166. "occurs.",
  167. -1: "N must be positive.",
  168. -2: "NEV must be positive.",
  169. -3: "NCV-NEV >= 1 and less than or equal to N.",
  170. -5: "WHICH must be one of 'LM', 'SM', 'LR', 'SR', 'LI', 'SI'",
  171. -6: "BMAT must be one of 'I' or 'G'.",
  172. -7: "Length of private work WORKL array is not sufficient.",
  173. -8: "Error return from LAPACK eigenvalue calculation. "
  174. "This should never happened.",
  175. -9: "Error return from calculation of eigenvectors. "
  176. "Informational error from LAPACK routine ztrevc.",
  177. -10: "IPARAM(7) must be 1,2,3",
  178. -11: "IPARAM(7) = 1 and BMAT = 'G' are incompatible.",
  179. -12: "HOWMNY = 'S' not yet implemented",
  180. -13: "HOWMNY must be one of 'A' or 'P' if RVEC = .true.",
  181. -14: "ZNAUPD did not find any eigenvalues to sufficient "
  182. "accuracy.",
  183. -15: "ZNEUPD got a different count of the number of "
  184. "converged Ritz values than ZNAUPD got. This "
  185. "indicates the user probably made an error in passing "
  186. "data from ZNAUPD to ZNEUPD or that the data was "
  187. "modified before entering ZNEUPD"
  188. }
  189. CNEUPD_ERRORS = ZNEUPD_ERRORS.copy()
  190. CNEUPD_ERRORS[-14] = ("CNAUPD did not find any eigenvalues to sufficient "
  191. "accuracy.")
  192. CNEUPD_ERRORS[-15] = ("CNEUPD got a different count of the number of "
  193. "converged Ritz values than CNAUPD got. This indicates "
  194. "the user probably made an error in passing data from "
  195. "CNAUPD to CNEUPD or that the data was modified before "
  196. "entering CNEUPD")
  197. DSEUPD_ERRORS = {
  198. 0: "Normal exit.",
  199. -1: "N must be positive.",
  200. -2: "NEV must be positive.",
  201. -3: "NCV must be greater than NEV and less than or equal to N.",
  202. -5: "WHICH must be one of 'LM', 'SM', 'LA', 'SA' or 'BE'.",
  203. -6: "BMAT must be one of 'I' or 'G'.",
  204. -7: "Length of private work WORKL array is not sufficient.",
  205. -8: ("Error return from trid. eigenvalue calculation; "
  206. "Information error from LAPACK routine dsteqr."),
  207. -9: "Starting vector is zero.",
  208. -10: "IPARAM(7) must be 1,2,3,4,5.",
  209. -11: "IPARAM(7) = 1 and BMAT = 'G' are incompatible.",
  210. -12: "NEV and WHICH = 'BE' are incompatible.",
  211. -14: "DSAUPD did not find any eigenvalues to sufficient accuracy.",
  212. -15: "HOWMNY must be one of 'A' or 'S' if RVEC = .true.",
  213. -16: "HOWMNY = 'S' not yet implemented",
  214. -17: ("DSEUPD got a different count of the number of converged "
  215. "Ritz values than DSAUPD got. This indicates the user "
  216. "probably made an error in passing data from DSAUPD to "
  217. "DSEUPD or that the data was modified before entering "
  218. "DSEUPD.")
  219. }
  220. SSEUPD_ERRORS = DSEUPD_ERRORS.copy()
  221. SSEUPD_ERRORS[-14] = ("SSAUPD did not find any eigenvalues "
  222. "to sufficient accuracy.")
  223. SSEUPD_ERRORS[-17] = ("SSEUPD got a different count of the number of "
  224. "converged "
  225. "Ritz values than SSAUPD got. This indicates the user "
  226. "probably made an error in passing data from SSAUPD to "
  227. "SSEUPD or that the data was modified before entering "
  228. "SSEUPD.")
  229. _SAUPD_ERRORS = {'d': DSAUPD_ERRORS,
  230. 's': SSAUPD_ERRORS}
  231. _NAUPD_ERRORS = {'d': DNAUPD_ERRORS,
  232. 's': SNAUPD_ERRORS,
  233. 'z': ZNAUPD_ERRORS,
  234. 'c': CNAUPD_ERRORS}
  235. _SEUPD_ERRORS = {'d': DSEUPD_ERRORS,
  236. 's': SSEUPD_ERRORS}
  237. _NEUPD_ERRORS = {'d': DNEUPD_ERRORS,
  238. 's': SNEUPD_ERRORS,
  239. 'z': ZNEUPD_ERRORS,
  240. 'c': CNEUPD_ERRORS}
  241. # accepted values of parameter WHICH in _SEUPD
  242. _SEUPD_WHICH = ['LM', 'SM', 'LA', 'SA', 'BE']
  243. # accepted values of parameter WHICH in _NAUPD
  244. _NEUPD_WHICH = ['LM', 'SM', 'LR', 'SR', 'LI', 'SI']
  245. class ArpackError(RuntimeError):
  246. """
  247. ARPACK error
  248. """
  249. def __init__(self, info, infodict=_NAUPD_ERRORS):
  250. msg = infodict.get(info, "Unknown error")
  251. RuntimeError.__init__(self, "ARPACK error %d: %s" % (info, msg))
  252. class ArpackNoConvergence(ArpackError):
  253. """
  254. ARPACK iteration did not converge
  255. Attributes
  256. ----------
  257. eigenvalues : ndarray
  258. Partial result. Converged eigenvalues.
  259. eigenvectors : ndarray
  260. Partial result. Converged eigenvectors.
  261. """
  262. def __init__(self, msg, eigenvalues, eigenvectors):
  263. ArpackError.__init__(self, -1, {-1: msg})
  264. self.eigenvalues = eigenvalues
  265. self.eigenvectors = eigenvectors
  266. def choose_ncv(k):
  267. """
  268. Choose number of lanczos vectors based on target number
  269. of singular/eigen values and vectors to compute, k.
  270. """
  271. return max(2 * k + 1, 20)
  272. class _ArpackParams(object):
  273. def __init__(self, n, k, tp, mode=1, sigma=None,
  274. ncv=None, v0=None, maxiter=None, which="LM", tol=0):
  275. if k <= 0:
  276. raise ValueError("k must be positive, k=%d" % k)
  277. if maxiter is None:
  278. maxiter = n * 10
  279. if maxiter <= 0:
  280. raise ValueError("maxiter must be positive, maxiter=%d" % maxiter)
  281. if tp not in 'fdFD':
  282. raise ValueError("matrix type must be 'f', 'd', 'F', or 'D'")
  283. if v0 is not None:
  284. # ARPACK overwrites its initial resid, make a copy
  285. self.resid = np.array(v0, copy=True)
  286. info = 1
  287. else:
  288. # ARPACK will use a random initial vector.
  289. self.resid = np.zeros(n, tp)
  290. info = 0
  291. if sigma is None:
  292. #sigma not used
  293. self.sigma = 0
  294. else:
  295. self.sigma = sigma
  296. if ncv is None:
  297. ncv = choose_ncv(k)
  298. ncv = min(ncv, n)
  299. self.v = np.zeros((n, ncv), tp) # holds Ritz vectors
  300. self.iparam = np.zeros(11, "int")
  301. # set solver mode and parameters
  302. ishfts = 1
  303. self.mode = mode
  304. self.iparam[0] = ishfts
  305. self.iparam[2] = maxiter
  306. self.iparam[3] = 1
  307. self.iparam[6] = mode
  308. self.n = n
  309. self.tol = tol
  310. self.k = k
  311. self.maxiter = maxiter
  312. self.ncv = ncv
  313. self.which = which
  314. self.tp = tp
  315. self.info = info
  316. self.converged = False
  317. self.ido = 0
  318. def _raise_no_convergence(self):
  319. msg = "No convergence (%d iterations, %d/%d eigenvectors converged)"
  320. k_ok = self.iparam[4]
  321. num_iter = self.iparam[2]
  322. try:
  323. ev, vec = self.extract(True)
  324. except ArpackError as err:
  325. msg = "%s [%s]" % (msg, err)
  326. ev = np.zeros((0,))
  327. vec = np.zeros((self.n, 0))
  328. k_ok = 0
  329. raise ArpackNoConvergence(msg % (num_iter, k_ok, self.k), ev, vec)
  330. class _SymmetricArpackParams(_ArpackParams):
  331. def __init__(self, n, k, tp, matvec, mode=1, M_matvec=None,
  332. Minv_matvec=None, sigma=None,
  333. ncv=None, v0=None, maxiter=None, which="LM", tol=0):
  334. # The following modes are supported:
  335. # mode = 1:
  336. # Solve the standard eigenvalue problem:
  337. # A*x = lambda*x :
  338. # A - symmetric
  339. # Arguments should be
  340. # matvec = left multiplication by A
  341. # M_matvec = None [not used]
  342. # Minv_matvec = None [not used]
  343. #
  344. # mode = 2:
  345. # Solve the general eigenvalue problem:
  346. # A*x = lambda*M*x
  347. # A - symmetric
  348. # M - symmetric positive definite
  349. # Arguments should be
  350. # matvec = left multiplication by A
  351. # M_matvec = left multiplication by M
  352. # Minv_matvec = left multiplication by M^-1
  353. #
  354. # mode = 3:
  355. # Solve the general eigenvalue problem in shift-invert mode:
  356. # A*x = lambda*M*x
  357. # A - symmetric
  358. # M - symmetric positive semi-definite
  359. # Arguments should be
  360. # matvec = None [not used]
  361. # M_matvec = left multiplication by M
  362. # or None, if M is the identity
  363. # Minv_matvec = left multiplication by [A-sigma*M]^-1
  364. #
  365. # mode = 4:
  366. # Solve the general eigenvalue problem in Buckling mode:
  367. # A*x = lambda*AG*x
  368. # A - symmetric positive semi-definite
  369. # AG - symmetric indefinite
  370. # Arguments should be
  371. # matvec = left multiplication by A
  372. # M_matvec = None [not used]
  373. # Minv_matvec = left multiplication by [A-sigma*AG]^-1
  374. #
  375. # mode = 5:
  376. # Solve the general eigenvalue problem in Cayley-transformed mode:
  377. # A*x = lambda*M*x
  378. # A - symmetric
  379. # M - symmetric positive semi-definite
  380. # Arguments should be
  381. # matvec = left multiplication by A
  382. # M_matvec = left multiplication by M
  383. # or None, if M is the identity
  384. # Minv_matvec = left multiplication by [A-sigma*M]^-1
  385. if mode == 1:
  386. if matvec is None:
  387. raise ValueError("matvec must be specified for mode=1")
  388. if M_matvec is not None:
  389. raise ValueError("M_matvec cannot be specified for mode=1")
  390. if Minv_matvec is not None:
  391. raise ValueError("Minv_matvec cannot be specified for mode=1")
  392. self.OP = matvec
  393. self.B = lambda x: x
  394. self.bmat = 'I'
  395. elif mode == 2:
  396. if matvec is None:
  397. raise ValueError("matvec must be specified for mode=2")
  398. if M_matvec is None:
  399. raise ValueError("M_matvec must be specified for mode=2")
  400. if Minv_matvec is None:
  401. raise ValueError("Minv_matvec must be specified for mode=2")
  402. self.OP = lambda x: Minv_matvec(matvec(x))
  403. self.OPa = Minv_matvec
  404. self.OPb = matvec
  405. self.B = M_matvec
  406. self.bmat = 'G'
  407. elif mode == 3:
  408. if matvec is not None:
  409. raise ValueError("matvec must not be specified for mode=3")
  410. if Minv_matvec is None:
  411. raise ValueError("Minv_matvec must be specified for mode=3")
  412. if M_matvec is None:
  413. self.OP = Minv_matvec
  414. self.OPa = Minv_matvec
  415. self.B = lambda x: x
  416. self.bmat = 'I'
  417. else:
  418. self.OP = lambda x: Minv_matvec(M_matvec(x))
  419. self.OPa = Minv_matvec
  420. self.B = M_matvec
  421. self.bmat = 'G'
  422. elif mode == 4:
  423. if matvec is None:
  424. raise ValueError("matvec must be specified for mode=4")
  425. if M_matvec is not None:
  426. raise ValueError("M_matvec must not be specified for mode=4")
  427. if Minv_matvec is None:
  428. raise ValueError("Minv_matvec must be specified for mode=4")
  429. self.OPa = Minv_matvec
  430. self.OP = lambda x: self.OPa(matvec(x))
  431. self.B = matvec
  432. self.bmat = 'G'
  433. elif mode == 5:
  434. if matvec is None:
  435. raise ValueError("matvec must be specified for mode=5")
  436. if Minv_matvec is None:
  437. raise ValueError("Minv_matvec must be specified for mode=5")
  438. self.OPa = Minv_matvec
  439. self.A_matvec = matvec
  440. if M_matvec is None:
  441. self.OP = lambda x: Minv_matvec(matvec(x) + sigma * x)
  442. self.B = lambda x: x
  443. self.bmat = 'I'
  444. else:
  445. self.OP = lambda x: Minv_matvec(matvec(x)
  446. + sigma * M_matvec(x))
  447. self.B = M_matvec
  448. self.bmat = 'G'
  449. else:
  450. raise ValueError("mode=%i not implemented" % mode)
  451. if which not in _SEUPD_WHICH:
  452. raise ValueError("which must be one of %s"
  453. % ' '.join(_SEUPD_WHICH))
  454. if k >= n:
  455. raise ValueError("k must be less than ndim(A), k=%d" % k)
  456. _ArpackParams.__init__(self, n, k, tp, mode, sigma,
  457. ncv, v0, maxiter, which, tol)
  458. if self.ncv > n or self.ncv <= k:
  459. raise ValueError("ncv must be k<ncv<=n, ncv=%s" % self.ncv)
  460. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  461. self.workd = _aligned_zeros(3 * n, self.tp)
  462. self.workl = _aligned_zeros(self.ncv * (self.ncv + 8), self.tp)
  463. ltr = _type_conv[self.tp]
  464. if ltr not in ["s", "d"]:
  465. raise ValueError("Input matrix is not real-valued.")
  466. self._arpack_solver = _arpack.__dict__[ltr + 'saupd']
  467. self._arpack_extract = _arpack.__dict__[ltr + 'seupd']
  468. self.iterate_infodict = _SAUPD_ERRORS[ltr]
  469. self.extract_infodict = _SEUPD_ERRORS[ltr]
  470. self.ipntr = np.zeros(11, "int")
  471. def iterate(self):
  472. self.ido, self.tol, self.resid, self.v, self.iparam, self.ipntr, self.info = \
  473. self._arpack_solver(self.ido, self.bmat, self.which, self.k,
  474. self.tol, self.resid, self.v, self.iparam,
  475. self.ipntr, self.workd, self.workl, self.info)
  476. xslice = slice(self.ipntr[0] - 1, self.ipntr[0] - 1 + self.n)
  477. yslice = slice(self.ipntr[1] - 1, self.ipntr[1] - 1 + self.n)
  478. if self.ido == -1:
  479. # initialization
  480. self.workd[yslice] = self.OP(self.workd[xslice])
  481. elif self.ido == 1:
  482. # compute y = Op*x
  483. if self.mode == 1:
  484. self.workd[yslice] = self.OP(self.workd[xslice])
  485. elif self.mode == 2:
  486. self.workd[xslice] = self.OPb(self.workd[xslice])
  487. self.workd[yslice] = self.OPa(self.workd[xslice])
  488. elif self.mode == 5:
  489. Bxslice = slice(self.ipntr[2] - 1, self.ipntr[2] - 1 + self.n)
  490. Ax = self.A_matvec(self.workd[xslice])
  491. self.workd[yslice] = self.OPa(Ax + (self.sigma *
  492. self.workd[Bxslice]))
  493. else:
  494. Bxslice = slice(self.ipntr[2] - 1, self.ipntr[2] - 1 + self.n)
  495. self.workd[yslice] = self.OPa(self.workd[Bxslice])
  496. elif self.ido == 2:
  497. self.workd[yslice] = self.B(self.workd[xslice])
  498. elif self.ido == 3:
  499. raise ValueError("ARPACK requested user shifts. Assure ISHIFT==0")
  500. else:
  501. self.converged = True
  502. if self.info == 0:
  503. pass
  504. elif self.info == 1:
  505. self._raise_no_convergence()
  506. else:
  507. raise ArpackError(self.info, infodict=self.iterate_infodict)
  508. def extract(self, return_eigenvectors):
  509. rvec = return_eigenvectors
  510. ierr = 0
  511. howmny = 'A' # return all eigenvectors
  512. sselect = np.zeros(self.ncv, 'int') # unused
  513. d, z, ierr = self._arpack_extract(rvec, howmny, sselect, self.sigma,
  514. self.bmat, self.which, self.k,
  515. self.tol, self.resid, self.v,
  516. self.iparam[0:7], self.ipntr,
  517. self.workd[0:2 * self.n],
  518. self.workl, ierr)
  519. if ierr != 0:
  520. raise ArpackError(ierr, infodict=self.extract_infodict)
  521. k_ok = self.iparam[4]
  522. d = d[:k_ok]
  523. z = z[:, :k_ok]
  524. if return_eigenvectors:
  525. return d, z
  526. else:
  527. return d
  528. class _UnsymmetricArpackParams(_ArpackParams):
  529. def __init__(self, n, k, tp, matvec, mode=1, M_matvec=None,
  530. Minv_matvec=None, sigma=None,
  531. ncv=None, v0=None, maxiter=None, which="LM", tol=0):
  532. # The following modes are supported:
  533. # mode = 1:
  534. # Solve the standard eigenvalue problem:
  535. # A*x = lambda*x
  536. # A - square matrix
  537. # Arguments should be
  538. # matvec = left multiplication by A
  539. # M_matvec = None [not used]
  540. # Minv_matvec = None [not used]
  541. #
  542. # mode = 2:
  543. # Solve the generalized eigenvalue problem:
  544. # A*x = lambda*M*x
  545. # A - square matrix
  546. # M - symmetric, positive semi-definite
  547. # Arguments should be
  548. # matvec = left multiplication by A
  549. # M_matvec = left multiplication by M
  550. # Minv_matvec = left multiplication by M^-1
  551. #
  552. # mode = 3,4:
  553. # Solve the general eigenvalue problem in shift-invert mode:
  554. # A*x = lambda*M*x
  555. # A - square matrix
  556. # M - symmetric, positive semi-definite
  557. # Arguments should be
  558. # matvec = None [not used]
  559. # M_matvec = left multiplication by M
  560. # or None, if M is the identity
  561. # Minv_matvec = left multiplication by [A-sigma*M]^-1
  562. # if A is real and mode==3, use the real part of Minv_matvec
  563. # if A is real and mode==4, use the imag part of Minv_matvec
  564. # if A is complex and mode==3,
  565. # use real and imag parts of Minv_matvec
  566. if mode == 1:
  567. if matvec is None:
  568. raise ValueError("matvec must be specified for mode=1")
  569. if M_matvec is not None:
  570. raise ValueError("M_matvec cannot be specified for mode=1")
  571. if Minv_matvec is not None:
  572. raise ValueError("Minv_matvec cannot be specified for mode=1")
  573. self.OP = matvec
  574. self.B = lambda x: x
  575. self.bmat = 'I'
  576. elif mode == 2:
  577. if matvec is None:
  578. raise ValueError("matvec must be specified for mode=2")
  579. if M_matvec is None:
  580. raise ValueError("M_matvec must be specified for mode=2")
  581. if Minv_matvec is None:
  582. raise ValueError("Minv_matvec must be specified for mode=2")
  583. self.OP = lambda x: Minv_matvec(matvec(x))
  584. self.OPa = Minv_matvec
  585. self.OPb = matvec
  586. self.B = M_matvec
  587. self.bmat = 'G'
  588. elif mode in (3, 4):
  589. if matvec is None:
  590. raise ValueError("matvec must be specified "
  591. "for mode in (3,4)")
  592. if Minv_matvec is None:
  593. raise ValueError("Minv_matvec must be specified "
  594. "for mode in (3,4)")
  595. self.matvec = matvec
  596. if tp in 'DF': # complex type
  597. if mode == 3:
  598. self.OPa = Minv_matvec
  599. else:
  600. raise ValueError("mode=4 invalid for complex A")
  601. else: # real type
  602. if mode == 3:
  603. self.OPa = lambda x: np.real(Minv_matvec(x))
  604. else:
  605. self.OPa = lambda x: np.imag(Minv_matvec(x))
  606. if M_matvec is None:
  607. self.B = lambda x: x
  608. self.bmat = 'I'
  609. self.OP = self.OPa
  610. else:
  611. self.B = M_matvec
  612. self.bmat = 'G'
  613. self.OP = lambda x: self.OPa(M_matvec(x))
  614. else:
  615. raise ValueError("mode=%i not implemented" % mode)
  616. if which not in _NEUPD_WHICH:
  617. raise ValueError("Parameter which must be one of %s"
  618. % ' '.join(_NEUPD_WHICH))
  619. if k >= n - 1:
  620. raise ValueError("k must be less than ndim(A)-1, k=%d" % k)
  621. _ArpackParams.__init__(self, n, k, tp, mode, sigma,
  622. ncv, v0, maxiter, which, tol)
  623. if self.ncv > n or self.ncv <= k + 1:
  624. raise ValueError("ncv must be k+1<ncv<=n, ncv=%s" % self.ncv)
  625. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  626. self.workd = _aligned_zeros(3 * n, self.tp)
  627. self.workl = _aligned_zeros(3 * self.ncv * (self.ncv + 2), self.tp)
  628. ltr = _type_conv[self.tp]
  629. self._arpack_solver = _arpack.__dict__[ltr + 'naupd']
  630. self._arpack_extract = _arpack.__dict__[ltr + 'neupd']
  631. self.iterate_infodict = _NAUPD_ERRORS[ltr]
  632. self.extract_infodict = _NEUPD_ERRORS[ltr]
  633. self.ipntr = np.zeros(14, "int")
  634. if self.tp in 'FD':
  635. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  636. self.rwork = _aligned_zeros(self.ncv, self.tp.lower())
  637. else:
  638. self.rwork = None
  639. def iterate(self):
  640. if self.tp in 'fd':
  641. self.ido, self.tol, self.resid, self.v, self.iparam, self.ipntr, self.info =\
  642. self._arpack_solver(self.ido, self.bmat, self.which, self.k,
  643. self.tol, self.resid, self.v, self.iparam,
  644. self.ipntr, self.workd, self.workl,
  645. self.info)
  646. else:
  647. self.ido, self.tol, self.resid, self.v, self.iparam, self.ipntr, self.info =\
  648. self._arpack_solver(self.ido, self.bmat, self.which, self.k,
  649. self.tol, self.resid, self.v, self.iparam,
  650. self.ipntr, self.workd, self.workl,
  651. self.rwork, self.info)
  652. xslice = slice(self.ipntr[0] - 1, self.ipntr[0] - 1 + self.n)
  653. yslice = slice(self.ipntr[1] - 1, self.ipntr[1] - 1 + self.n)
  654. if self.ido == -1:
  655. # initialization
  656. self.workd[yslice] = self.OP(self.workd[xslice])
  657. elif self.ido == 1:
  658. # compute y = Op*x
  659. if self.mode in (1, 2):
  660. self.workd[yslice] = self.OP(self.workd[xslice])
  661. else:
  662. Bxslice = slice(self.ipntr[2] - 1, self.ipntr[2] - 1 + self.n)
  663. self.workd[yslice] = self.OPa(self.workd[Bxslice])
  664. elif self.ido == 2:
  665. self.workd[yslice] = self.B(self.workd[xslice])
  666. elif self.ido == 3:
  667. raise ValueError("ARPACK requested user shifts. Assure ISHIFT==0")
  668. else:
  669. self.converged = True
  670. if self.info == 0:
  671. pass
  672. elif self.info == 1:
  673. self._raise_no_convergence()
  674. else:
  675. raise ArpackError(self.info, infodict=self.iterate_infodict)
  676. def extract(self, return_eigenvectors):
  677. k, n = self.k, self.n
  678. ierr = 0
  679. howmny = 'A' # return all eigenvectors
  680. sselect = np.zeros(self.ncv, 'int') # unused
  681. sigmar = np.real(self.sigma)
  682. sigmai = np.imag(self.sigma)
  683. workev = np.zeros(3 * self.ncv, self.tp)
  684. if self.tp in 'fd':
  685. dr = np.zeros(k + 1, self.tp)
  686. di = np.zeros(k + 1, self.tp)
  687. zr = np.zeros((n, k + 1), self.tp)
  688. dr, di, zr, ierr = \
  689. self._arpack_extract(return_eigenvectors,
  690. howmny, sselect, sigmar, sigmai, workev,
  691. self.bmat, self.which, k, self.tol, self.resid,
  692. self.v, self.iparam, self.ipntr,
  693. self.workd, self.workl, self.info)
  694. if ierr != 0:
  695. raise ArpackError(ierr, infodict=self.extract_infodict)
  696. nreturned = self.iparam[4] # number of good eigenvalues returned
  697. # Build complex eigenvalues from real and imaginary parts
  698. d = dr + 1.0j * di
  699. # Arrange the eigenvectors: complex eigenvectors are stored as
  700. # real,imaginary in consecutive columns
  701. z = zr.astype(self.tp.upper())
  702. # The ARPACK nonsymmetric real and double interface (s,d)naupd
  703. # return eigenvalues and eigenvectors in real (float,double)
  704. # arrays.
  705. # Efficiency: this should check that return_eigenvectors == True
  706. # before going through this construction.
  707. if sigmai == 0:
  708. i = 0
  709. while i <= k:
  710. # check if complex
  711. if abs(d[i].imag) != 0:
  712. # this is a complex conjugate pair with eigenvalues
  713. # in consecutive columns
  714. if i < k:
  715. z[:, i] = zr[:, i] + 1.0j * zr[:, i + 1]
  716. z[:, i + 1] = z[:, i].conjugate()
  717. i += 1
  718. else:
  719. #last eigenvalue is complex: the imaginary part of
  720. # the eigenvector has not been returned
  721. #this can only happen if nreturned > k, so we'll
  722. # throw out this case.
  723. nreturned -= 1
  724. i += 1
  725. else:
  726. # real matrix, mode 3 or 4, imag(sigma) is nonzero:
  727. # see remark 3 in <s,d>neupd.f
  728. # Build complex eigenvalues from real and imaginary parts
  729. i = 0
  730. while i <= k:
  731. if abs(d[i].imag) == 0:
  732. d[i] = np.dot(zr[:, i], self.matvec(zr[:, i]))
  733. else:
  734. if i < k:
  735. z[:, i] = zr[:, i] + 1.0j * zr[:, i + 1]
  736. z[:, i + 1] = z[:, i].conjugate()
  737. d[i] = ((np.dot(zr[:, i],
  738. self.matvec(zr[:, i]))
  739. + np.dot(zr[:, i + 1],
  740. self.matvec(zr[:, i + 1])))
  741. + 1j * (np.dot(zr[:, i],
  742. self.matvec(zr[:, i + 1]))
  743. - np.dot(zr[:, i + 1],
  744. self.matvec(zr[:, i]))))
  745. d[i + 1] = d[i].conj()
  746. i += 1
  747. else:
  748. #last eigenvalue is complex: the imaginary part of
  749. # the eigenvector has not been returned
  750. #this can only happen if nreturned > k, so we'll
  751. # throw out this case.
  752. nreturned -= 1
  753. i += 1
  754. # Now we have k+1 possible eigenvalues and eigenvectors
  755. # Return the ones specified by the keyword "which"
  756. if nreturned <= k:
  757. # we got less or equal as many eigenvalues we wanted
  758. d = d[:nreturned]
  759. z = z[:, :nreturned]
  760. else:
  761. # we got one extra eigenvalue (likely a cc pair, but which?)
  762. if self.mode in (1, 2):
  763. rd = d
  764. elif self.mode in (3, 4):
  765. rd = 1 / (d - self.sigma)
  766. if self.which in ['LR', 'SR']:
  767. ind = np.argsort(rd.real)
  768. elif self.which in ['LI', 'SI']:
  769. # for LI,SI ARPACK returns largest,smallest
  770. # abs(imaginary) (complex pairs come together)
  771. ind = np.argsort(abs(rd.imag))
  772. else:
  773. ind = np.argsort(abs(rd))
  774. if self.which in ['LR', 'LM', 'LI']:
  775. ind = ind[-k:][::-1]
  776. elif self.which in ['SR', 'SM', 'SI']:
  777. ind = ind[:k]
  778. d = d[ind]
  779. z = z[:, ind]
  780. else:
  781. # complex is so much simpler...
  782. d, z, ierr =\
  783. self._arpack_extract(return_eigenvectors,
  784. howmny, sselect, self.sigma, workev,
  785. self.bmat, self.which, k, self.tol, self.resid,
  786. self.v, self.iparam, self.ipntr,
  787. self.workd, self.workl, self.rwork, ierr)
  788. if ierr != 0:
  789. raise ArpackError(ierr, infodict=self.extract_infodict)
  790. k_ok = self.iparam[4]
  791. d = d[:k_ok]
  792. z = z[:, :k_ok]
  793. if return_eigenvectors:
  794. return d, z
  795. else:
  796. return d
  797. def _aslinearoperator_with_dtype(m):
  798. m = aslinearoperator(m)
  799. if not hasattr(m, 'dtype'):
  800. x = np.zeros(m.shape[1])
  801. m.dtype = (m * x).dtype
  802. return m
  803. class SpLuInv(LinearOperator):
  804. """
  805. SpLuInv:
  806. helper class to repeatedly solve M*x=b
  807. using a sparse LU-decopposition of M
  808. """
  809. def __init__(self, M):
  810. self.M_lu = splu(M)
  811. self.shape = M.shape
  812. self.dtype = M.dtype
  813. self.isreal = not np.issubdtype(self.dtype, np.complexfloating)
  814. def _matvec(self, x):
  815. # careful here: splu.solve will throw away imaginary
  816. # part of x if M is real
  817. x = np.asarray(x)
  818. if self.isreal and np.issubdtype(x.dtype, np.complexfloating):
  819. return (self.M_lu.solve(np.real(x).astype(self.dtype))
  820. + 1j * self.M_lu.solve(np.imag(x).astype(self.dtype)))
  821. else:
  822. return self.M_lu.solve(x.astype(self.dtype))
  823. class LuInv(LinearOperator):
  824. """
  825. LuInv:
  826. helper class to repeatedly solve M*x=b
  827. using an LU-decomposition of M
  828. """
  829. def __init__(self, M):
  830. self.M_lu = lu_factor(M)
  831. self.shape = M.shape
  832. self.dtype = M.dtype
  833. def _matvec(self, x):
  834. return lu_solve(self.M_lu, x)
  835. def gmres_loose(A, b, tol):
  836. """
  837. gmres with looser termination condition.
  838. """
  839. b = np.asarray(b)
  840. min_tol = 1000 * np.sqrt(b.size) * np.finfo(b.dtype).eps
  841. return gmres(A, b, tol=max(tol, min_tol), atol=0)
  842. class IterInv(LinearOperator):
  843. """
  844. IterInv:
  845. helper class to repeatedly solve M*x=b
  846. using an iterative method.
  847. """
  848. def __init__(self, M, ifunc=gmres_loose, tol=0):
  849. self.M = M
  850. if hasattr(M, 'dtype'):
  851. self.dtype = M.dtype
  852. else:
  853. x = np.zeros(M.shape[1])
  854. self.dtype = (M * x).dtype
  855. self.shape = M.shape
  856. if tol <= 0:
  857. # when tol=0, ARPACK uses machine tolerance as calculated
  858. # by LAPACK's _LAMCH function. We should match this
  859. tol = 2 * np.finfo(self.dtype).eps
  860. self.ifunc = ifunc
  861. self.tol = tol
  862. def _matvec(self, x):
  863. b, info = self.ifunc(self.M, x, tol=self.tol)
  864. if info != 0:
  865. raise ValueError("Error in inverting M: function "
  866. "%s did not converge (info = %i)."
  867. % (self.ifunc.__name__, info))
  868. return b
  869. class IterOpInv(LinearOperator):
  870. """
  871. IterOpInv:
  872. helper class to repeatedly solve [A-sigma*M]*x = b
  873. using an iterative method
  874. """
  875. def __init__(self, A, M, sigma, ifunc=gmres_loose, tol=0):
  876. self.A = A
  877. self.M = M
  878. self.sigma = sigma
  879. def mult_func(x):
  880. return A.matvec(x) - sigma * M.matvec(x)
  881. def mult_func_M_None(x):
  882. return A.matvec(x) - sigma * x
  883. x = np.zeros(A.shape[1])
  884. if M is None:
  885. dtype = mult_func_M_None(x).dtype
  886. self.OP = LinearOperator(self.A.shape,
  887. mult_func_M_None,
  888. dtype=dtype)
  889. else:
  890. dtype = mult_func(x).dtype
  891. self.OP = LinearOperator(self.A.shape,
  892. mult_func,
  893. dtype=dtype)
  894. self.shape = A.shape
  895. if tol <= 0:
  896. # when tol=0, ARPACK uses machine tolerance as calculated
  897. # by LAPACK's _LAMCH function. We should match this
  898. tol = 2 * np.finfo(self.OP.dtype).eps
  899. self.ifunc = ifunc
  900. self.tol = tol
  901. def _matvec(self, x):
  902. b, info = self.ifunc(self.OP, x, tol=self.tol)
  903. if info != 0:
  904. raise ValueError("Error in inverting [A-sigma*M]: function "
  905. "%s did not converge (info = %i)."
  906. % (self.ifunc.__name__, info))
  907. return b
  908. @property
  909. def dtype(self):
  910. return self.OP.dtype
  911. def _fast_spmatrix_to_csc(A, hermitian=False):
  912. """Convert sparse matrix to CSC (by transposing, if possible)"""
  913. if (isspmatrix_csr(A) and hermitian
  914. and not np.issubdtype(A.dtype, np.complexfloating)):
  915. return A.T
  916. else:
  917. return A.tocsc()
  918. def get_inv_matvec(M, hermitian=False, tol=0):
  919. if isdense(M):
  920. return LuInv(M).matvec
  921. elif isspmatrix(M):
  922. M = _fast_spmatrix_to_csc(M, hermitian=hermitian)
  923. return SpLuInv(M).matvec
  924. else:
  925. return IterInv(M, tol=tol).matvec
  926. def get_OPinv_matvec(A, M, sigma, hermitian=False, tol=0):
  927. if sigma == 0:
  928. return get_inv_matvec(A, hermitian=hermitian, tol=tol)
  929. if M is None:
  930. #M is the identity matrix
  931. if isdense(A):
  932. if (np.issubdtype(A.dtype, np.complexfloating)
  933. or np.imag(sigma) == 0):
  934. A = np.copy(A)
  935. else:
  936. A = A + 0j
  937. A.flat[::A.shape[1] + 1] -= sigma
  938. return LuInv(A).matvec
  939. elif isspmatrix(A):
  940. A = A - sigma * eye(A.shape[0])
  941. A = _fast_spmatrix_to_csc(A, hermitian=hermitian)
  942. return SpLuInv(A).matvec
  943. else:
  944. return IterOpInv(_aslinearoperator_with_dtype(A),
  945. M, sigma, tol=tol).matvec
  946. else:
  947. if ((not isdense(A) and not isspmatrix(A)) or
  948. (not isdense(M) and not isspmatrix(M))):
  949. return IterOpInv(_aslinearoperator_with_dtype(A),
  950. _aslinearoperator_with_dtype(M),
  951. sigma, tol=tol).matvec
  952. elif isdense(A) or isdense(M):
  953. return LuInv(A - sigma * M).matvec
  954. else:
  955. OP = A - sigma * M
  956. OP = _fast_spmatrix_to_csc(OP, hermitian=hermitian)
  957. return SpLuInv(OP).matvec
  958. # ARPACK is not threadsafe or reentrant (SAVE variables), so we need a
  959. # lock and a re-entering check.
  960. _ARPACK_LOCK = ReentrancyLock("Nested calls to eigs/eighs not allowed: "
  961. "ARPACK is not re-entrant")
  962. def eigs(A, k=6, M=None, sigma=None, which='LM', v0=None,
  963. ncv=None, maxiter=None, tol=0, return_eigenvectors=True,
  964. Minv=None, OPinv=None, OPpart=None):
  965. """
  966. Find k eigenvalues and eigenvectors of the square matrix A.
  967. Solves ``A * x[i] = w[i] * x[i]``, the standard eigenvalue problem
  968. for w[i] eigenvalues with corresponding eigenvectors x[i].
  969. If M is specified, solves ``A * x[i] = w[i] * M * x[i]``, the
  970. generalized eigenvalue problem for w[i] eigenvalues
  971. with corresponding eigenvectors x[i]
  972. Parameters
  973. ----------
  974. A : ndarray, sparse matrix or LinearOperator
  975. An array, sparse matrix, or LinearOperator representing
  976. the operation ``A * x``, where A is a real or complex square matrix.
  977. k : int, optional
  978. The number of eigenvalues and eigenvectors desired.
  979. `k` must be smaller than N-1. It is not possible to compute all
  980. eigenvectors of a matrix.
  981. M : ndarray, sparse matrix or LinearOperator, optional
  982. An array, sparse matrix, or LinearOperator representing
  983. the operation M*x for the generalized eigenvalue problem
  984. A * x = w * M * x.
  985. M must represent a real, symmetric matrix if A is real, and must
  986. represent a complex, hermitian matrix if A is complex. For best
  987. results, the data type of M should be the same as that of A.
  988. Additionally:
  989. If `sigma` is None, M is positive definite
  990. If sigma is specified, M is positive semi-definite
  991. If sigma is None, eigs requires an operator to compute the solution
  992. of the linear equation ``M * x = b``. This is done internally via a
  993. (sparse) LU decomposition for an explicit matrix M, or via an
  994. iterative solver for a general linear operator. Alternatively,
  995. the user can supply the matrix or operator Minv, which gives
  996. ``x = Minv * b = M^-1 * b``.
  997. sigma : real or complex, optional
  998. Find eigenvalues near sigma using shift-invert mode. This requires
  999. an operator to compute the solution of the linear system
  1000. ``[A - sigma * M] * x = b``, where M is the identity matrix if
  1001. unspecified. This is computed internally via a (sparse) LU
  1002. decomposition for explicit matrices A & M, or via an iterative
  1003. solver if either A or M is a general linear operator.
  1004. Alternatively, the user can supply the matrix or operator OPinv,
  1005. which gives ``x = OPinv * b = [A - sigma * M]^-1 * b``.
  1006. For a real matrix A, shift-invert can either be done in imaginary
  1007. mode or real mode, specified by the parameter OPpart ('r' or 'i').
  1008. Note that when sigma is specified, the keyword 'which' (below)
  1009. refers to the shifted eigenvalues ``w'[i]`` where:
  1010. If A is real and OPpart == 'r' (default),
  1011. ``w'[i] = 1/2 * [1/(w[i]-sigma) + 1/(w[i]-conj(sigma))]``.
  1012. If A is real and OPpart == 'i',
  1013. ``w'[i] = 1/2i * [1/(w[i]-sigma) - 1/(w[i]-conj(sigma))]``.
  1014. If A is complex, ``w'[i] = 1/(w[i]-sigma)``.
  1015. v0 : ndarray, optional
  1016. Starting vector for iteration.
  1017. Default: random
  1018. ncv : int, optional
  1019. The number of Lanczos vectors generated
  1020. `ncv` must be greater than `k`; it is recommended that ``ncv > 2*k``.
  1021. Default: ``min(n, max(2*k + 1, 20))``
  1022. which : str, ['LM' | 'SM' | 'LR' | 'SR' | 'LI' | 'SI'], optional
  1023. Which `k` eigenvectors and eigenvalues to find:
  1024. 'LM' : largest magnitude
  1025. 'SM' : smallest magnitude
  1026. 'LR' : largest real part
  1027. 'SR' : smallest real part
  1028. 'LI' : largest imaginary part
  1029. 'SI' : smallest imaginary part
  1030. When sigma != None, 'which' refers to the shifted eigenvalues w'[i]
  1031. (see discussion in 'sigma', above). ARPACK is generally better
  1032. at finding large values than small values. If small eigenvalues are
  1033. desired, consider using shift-invert mode for better performance.
  1034. maxiter : int, optional
  1035. Maximum number of Arnoldi update iterations allowed
  1036. Default: ``n*10``
  1037. tol : float, optional
  1038. Relative accuracy for eigenvalues (stopping criterion)
  1039. The default value of 0 implies machine precision.
  1040. return_eigenvectors : bool, optional
  1041. Return eigenvectors (True) in addition to eigenvalues
  1042. Minv : ndarray, sparse matrix or LinearOperator, optional
  1043. See notes in M, above.
  1044. OPinv : ndarray, sparse matrix or LinearOperator, optional
  1045. See notes in sigma, above.
  1046. OPpart : {'r' or 'i'}, optional
  1047. See notes in sigma, above
  1048. Returns
  1049. -------
  1050. w : ndarray
  1051. Array of k eigenvalues.
  1052. v : ndarray
  1053. An array of `k` eigenvectors.
  1054. ``v[:, i]`` is the eigenvector corresponding to the eigenvalue w[i].
  1055. Raises
  1056. ------
  1057. ArpackNoConvergence
  1058. When the requested convergence is not obtained.
  1059. The currently converged eigenvalues and eigenvectors can be found
  1060. as ``eigenvalues`` and ``eigenvectors`` attributes of the exception
  1061. object.
  1062. See Also
  1063. --------
  1064. eigsh : eigenvalues and eigenvectors for symmetric matrix A
  1065. svds : singular value decomposition for a matrix A
  1066. Notes
  1067. -----
  1068. This function is a wrapper to the ARPACK [1]_ SNEUPD, DNEUPD, CNEUPD,
  1069. ZNEUPD, functions which use the Implicitly Restarted Arnoldi Method to
  1070. find the eigenvalues and eigenvectors [2]_.
  1071. References
  1072. ----------
  1073. .. [1] ARPACK Software, http://www.caam.rice.edu/software/ARPACK/
  1074. .. [2] R. B. Lehoucq, D. C. Sorensen, and C. Yang, ARPACK USERS GUIDE:
  1075. Solution of Large Scale Eigenvalue Problems by Implicitly Restarted
  1076. Arnoldi Methods. SIAM, Philadelphia, PA, 1998.
  1077. Examples
  1078. --------
  1079. Find 6 eigenvectors of the identity matrix:
  1080. >>> from scipy.sparse.linalg import eigs
  1081. >>> id = np.eye(13)
  1082. >>> vals, vecs = eigs(id, k=6)
  1083. >>> vals
  1084. array([ 1.+0.j, 1.+0.j, 1.+0.j, 1.+0.j, 1.+0.j, 1.+0.j])
  1085. >>> vecs.shape
  1086. (13, 6)
  1087. """
  1088. if A.shape[0] != A.shape[1]:
  1089. raise ValueError('expected square matrix (shape=%s)' % (A.shape,))
  1090. if M is not None:
  1091. if M.shape != A.shape:
  1092. raise ValueError('wrong M dimensions %s, should be %s'
  1093. % (M.shape, A.shape))
  1094. if np.dtype(M.dtype).char.lower() != np.dtype(A.dtype).char.lower():
  1095. warnings.warn('M does not have the same type precision as A. '
  1096. 'This may adversely affect ARPACK convergence')
  1097. n = A.shape[0]
  1098. if k <= 0:
  1099. raise ValueError("k=%d must be greater than 0." % k)
  1100. if k >= n - 1:
  1101. warnings.warn("k >= N - 1 for N * N square matrix. "
  1102. "Attempting to use scipy.linalg.eig instead.",
  1103. RuntimeWarning)
  1104. if issparse(A):
  1105. raise TypeError("Cannot use scipy.linalg.eig for sparse A with "
  1106. "k >= N - 1. Use scipy.linalg.eig(A.toarray()) or"
  1107. " reduce k.")
  1108. if isinstance(A, LinearOperator):
  1109. raise TypeError("Cannot use scipy.linalg.eig for LinearOperator "
  1110. "A with k >= N - 1.")
  1111. if isinstance(M, LinearOperator):
  1112. raise TypeError("Cannot use scipy.linalg.eig for LinearOperator "
  1113. "M with k >= N - 1.")
  1114. return eig(A, b=M, right=return_eigenvectors)
  1115. if sigma is None:
  1116. matvec = _aslinearoperator_with_dtype(A).matvec
  1117. if OPinv is not None:
  1118. raise ValueError("OPinv should not be specified "
  1119. "with sigma = None.")
  1120. if OPpart is not None:
  1121. raise ValueError("OPpart should not be specified with "
  1122. "sigma = None or complex A")
  1123. if M is None:
  1124. #standard eigenvalue problem
  1125. mode = 1
  1126. M_matvec = None
  1127. Minv_matvec = None
  1128. if Minv is not None:
  1129. raise ValueError("Minv should not be "
  1130. "specified with M = None.")
  1131. else:
  1132. #general eigenvalue problem
  1133. mode = 2
  1134. if Minv is None:
  1135. Minv_matvec = get_inv_matvec(M, hermitian=True, tol=tol)
  1136. else:
  1137. Minv = _aslinearoperator_with_dtype(Minv)
  1138. Minv_matvec = Minv.matvec
  1139. M_matvec = _aslinearoperator_with_dtype(M).matvec
  1140. else:
  1141. #sigma is not None: shift-invert mode
  1142. if np.issubdtype(A.dtype, np.complexfloating):
  1143. if OPpart is not None:
  1144. raise ValueError("OPpart should not be specified "
  1145. "with sigma=None or complex A")
  1146. mode = 3
  1147. elif OPpart is None or OPpart.lower() == 'r':
  1148. mode = 3
  1149. elif OPpart.lower() == 'i':
  1150. if np.imag(sigma) == 0:
  1151. raise ValueError("OPpart cannot be 'i' if sigma is real")
  1152. mode = 4
  1153. else:
  1154. raise ValueError("OPpart must be one of ('r','i')")
  1155. matvec = _aslinearoperator_with_dtype(A).matvec
  1156. if Minv is not None:
  1157. raise ValueError("Minv should not be specified when sigma is")
  1158. if OPinv is None:
  1159. Minv_matvec = get_OPinv_matvec(A, M, sigma,
  1160. hermitian=False, tol=tol)
  1161. else:
  1162. OPinv = _aslinearoperator_with_dtype(OPinv)
  1163. Minv_matvec = OPinv.matvec
  1164. if M is None:
  1165. M_matvec = None
  1166. else:
  1167. M_matvec = _aslinearoperator_with_dtype(M).matvec
  1168. params = _UnsymmetricArpackParams(n, k, A.dtype.char, matvec, mode,
  1169. M_matvec, Minv_matvec, sigma,
  1170. ncv, v0, maxiter, which, tol)
  1171. with _ARPACK_LOCK:
  1172. while not params.converged:
  1173. params.iterate()
  1174. return params.extract(return_eigenvectors)
  1175. def eigsh(A, k=6, M=None, sigma=None, which='LM', v0=None,
  1176. ncv=None, maxiter=None, tol=0, return_eigenvectors=True,
  1177. Minv=None, OPinv=None, mode='normal'):
  1178. """
  1179. Find k eigenvalues and eigenvectors of the real symmetric square matrix
  1180. or complex hermitian matrix A.
  1181. Solves ``A * x[i] = w[i] * x[i]``, the standard eigenvalue problem for
  1182. w[i] eigenvalues with corresponding eigenvectors x[i].
  1183. If M is specified, solves ``A * x[i] = w[i] * M * x[i]``, the
  1184. generalized eigenvalue problem for w[i] eigenvalues
  1185. with corresponding eigenvectors x[i].
  1186. Parameters
  1187. ----------
  1188. A : ndarray, sparse matrix or LinearOperator
  1189. A square operator representing the operation ``A * x``, where ``A`` is
  1190. real symmetric or complex hermitian. For buckling mode (see below)
  1191. ``A`` must additionally be positive-definite.
  1192. k : int, optional
  1193. The number of eigenvalues and eigenvectors desired.
  1194. `k` must be smaller than N. It is not possible to compute all
  1195. eigenvectors of a matrix.
  1196. Returns
  1197. -------
  1198. w : array
  1199. Array of k eigenvalues.
  1200. v : array
  1201. An array representing the `k` eigenvectors. The column ``v[:, i]`` is
  1202. the eigenvector corresponding to the eigenvalue ``w[i]``.
  1203. Other Parameters
  1204. ----------------
  1205. M : An N x N matrix, array, sparse matrix, or linear operator representing
  1206. the operation ``M @ x`` for the generalized eigenvalue problem
  1207. A @ x = w * M @ x.
  1208. M must represent a real, symmetric matrix if A is real, and must
  1209. represent a complex, hermitian matrix if A is complex. For best
  1210. results, the data type of M should be the same as that of A.
  1211. Additionally:
  1212. If sigma is None, M is symmetric positive definite.
  1213. If sigma is specified, M is symmetric positive semi-definite.
  1214. In buckling mode, M is symmetric indefinite.
  1215. If sigma is None, eigsh requires an operator to compute the solution
  1216. of the linear equation ``M @ x = b``. This is done internally via a
  1217. (sparse) LU decomposition for an explicit matrix M, or via an
  1218. iterative solver for a general linear operator. Alternatively,
  1219. the user can supply the matrix or operator Minv, which gives
  1220. ``x = Minv @ b = M^-1 @ b``.
  1221. sigma : real
  1222. Find eigenvalues near sigma using shift-invert mode. This requires
  1223. an operator to compute the solution of the linear system
  1224. ``[A - sigma * M] x = b``, where M is the identity matrix if
  1225. unspecified. This is computed internally via a (sparse) LU
  1226. decomposition for explicit matrices A & M, or via an iterative
  1227. solver if either A or M is a general linear operator.
  1228. Alternatively, the user can supply the matrix or operator OPinv,
  1229. which gives ``x = OPinv @ b = [A - sigma * M]^-1 @ b``.
  1230. Note that when sigma is specified, the keyword 'which' refers to
  1231. the shifted eigenvalues ``w'[i]`` where:
  1232. if mode == 'normal', ``w'[i] = 1 / (w[i] - sigma)``.
  1233. if mode == 'cayley', ``w'[i] = (w[i] + sigma) / (w[i] - sigma)``.
  1234. if mode == 'buckling', ``w'[i] = w[i] / (w[i] - sigma)``.
  1235. (see further discussion in 'mode' below)
  1236. v0 : ndarray, optional
  1237. Starting vector for iteration.
  1238. Default: random
  1239. ncv : int, optional
  1240. The number of Lanczos vectors generated ncv must be greater than k and
  1241. smaller than n; it is recommended that ``ncv > 2*k``.
  1242. Default: ``min(n, max(2*k + 1, 20))``
  1243. which : str ['LM' | 'SM' | 'LA' | 'SA' | 'BE']
  1244. If A is a complex hermitian matrix, 'BE' is invalid.
  1245. Which `k` eigenvectors and eigenvalues to find:
  1246. 'LM' : Largest (in magnitude) eigenvalues.
  1247. 'SM' : Smallest (in magnitude) eigenvalues.
  1248. 'LA' : Largest (algebraic) eigenvalues.
  1249. 'SA' : Smallest (algebraic) eigenvalues.
  1250. 'BE' : Half (k/2) from each end of the spectrum.
  1251. When k is odd, return one more (k/2+1) from the high end.
  1252. When sigma != None, 'which' refers to the shifted eigenvalues ``w'[i]``
  1253. (see discussion in 'sigma', above). ARPACK is generally better
  1254. at finding large values than small values. If small eigenvalues are
  1255. desired, consider using shift-invert mode for better performance.
  1256. maxiter : int, optional
  1257. Maximum number of Arnoldi update iterations allowed.
  1258. Default: ``n*10``
  1259. tol : float
  1260. Relative accuracy for eigenvalues (stopping criterion).
  1261. The default value of 0 implies machine precision.
  1262. Minv : N x N matrix, array, sparse matrix, or LinearOperator
  1263. See notes in M, above.
  1264. OPinv : N x N matrix, array, sparse matrix, or LinearOperator
  1265. See notes in sigma, above.
  1266. return_eigenvectors : bool
  1267. Return eigenvectors (True) in addition to eigenvalues. This value determines
  1268. the order in which eigenvalues are sorted. The sort order is also dependent on the `which` variable.
  1269. For which = 'LM' or 'SA':
  1270. If `return_eigenvectors` is True, eigenvalues are sorted by algebraic value.
  1271. If `return_eigenvectors` is False, eigenvalues are sorted by absolute value.
  1272. For which = 'BE' or 'LA':
  1273. eigenvalues are always sorted by algebraic value.
  1274. For which = 'SM':
  1275. If `return_eigenvectors` is True, eigenvalues are sorted by algebraic value.
  1276. If `return_eigenvectors` is False, eigenvalues are sorted by decreasing absolute value.
  1277. mode : string ['normal' | 'buckling' | 'cayley']
  1278. Specify strategy to use for shift-invert mode. This argument applies
  1279. only for real-valued A and sigma != None. For shift-invert mode,
  1280. ARPACK internally solves the eigenvalue problem
  1281. ``OP * x'[i] = w'[i] * B * x'[i]``
  1282. and transforms the resulting Ritz vectors x'[i] and Ritz values w'[i]
  1283. into the desired eigenvectors and eigenvalues of the problem
  1284. ``A * x[i] = w[i] * M * x[i]``.
  1285. The modes are as follows:
  1286. 'normal' :
  1287. OP = [A - sigma * M]^-1 @ M,
  1288. B = M,
  1289. w'[i] = 1 / (w[i] - sigma)
  1290. 'buckling' :
  1291. OP = [A - sigma * M]^-1 @ A,
  1292. B = A,
  1293. w'[i] = w[i] / (w[i] - sigma)
  1294. 'cayley' :
  1295. OP = [A - sigma * M]^-1 @ [A + sigma * M],
  1296. B = M,
  1297. w'[i] = (w[i] + sigma) / (w[i] - sigma)
  1298. The choice of mode will affect which eigenvalues are selected by
  1299. the keyword 'which', and can also impact the stability of
  1300. convergence (see [2] for a discussion).
  1301. Raises
  1302. ------
  1303. ArpackNoConvergence
  1304. When the requested convergence is not obtained.
  1305. The currently converged eigenvalues and eigenvectors can be found
  1306. as ``eigenvalues`` and ``eigenvectors`` attributes of the exception
  1307. object.
  1308. See Also
  1309. --------
  1310. eigs : eigenvalues and eigenvectors for a general (nonsymmetric) matrix A
  1311. svds : singular value decomposition for a matrix A
  1312. Notes
  1313. -----
  1314. This function is a wrapper to the ARPACK [1]_ SSEUPD and DSEUPD
  1315. functions which use the Implicitly Restarted Lanczos Method to
  1316. find the eigenvalues and eigenvectors [2]_.
  1317. References
  1318. ----------
  1319. .. [1] ARPACK Software, http://www.caam.rice.edu/software/ARPACK/
  1320. .. [2] R. B. Lehoucq, D. C. Sorensen, and C. Yang, ARPACK USERS GUIDE:
  1321. Solution of Large Scale Eigenvalue Problems by Implicitly Restarted
  1322. Arnoldi Methods. SIAM, Philadelphia, PA, 1998.
  1323. Examples
  1324. --------
  1325. >>> from scipy.sparse.linalg import eigsh
  1326. >>> identity = np.eye(13)
  1327. >>> eigenvalues, eigenvectors = eigsh(identity, k=6)
  1328. >>> eigenvalues
  1329. array([1., 1., 1., 1., 1., 1.])
  1330. >>> eigenvectors.shape
  1331. (13, 6)
  1332. """
  1333. # complex hermitian matrices should be solved with eigs
  1334. if np.issubdtype(A.dtype, np.complexfloating):
  1335. if mode != 'normal':
  1336. raise ValueError("mode=%s cannot be used with "
  1337. "complex matrix A" % mode)
  1338. if which == 'BE':
  1339. raise ValueError("which='BE' cannot be used with complex matrix A")
  1340. elif which == 'LA':
  1341. which = 'LR'
  1342. elif which == 'SA':
  1343. which = 'SR'
  1344. ret = eigs(A, k, M=M, sigma=sigma, which=which, v0=v0,
  1345. ncv=ncv, maxiter=maxiter, tol=tol,
  1346. return_eigenvectors=return_eigenvectors, Minv=Minv,
  1347. OPinv=OPinv)
  1348. if return_eigenvectors:
  1349. return ret[0].real, ret[1]
  1350. else:
  1351. return ret.real
  1352. if A.shape[0] != A.shape[1]:
  1353. raise ValueError('expected square matrix (shape=%s)' % (A.shape,))
  1354. if M is not None:
  1355. if M.shape != A.shape:
  1356. raise ValueError('wrong M dimensions %s, should be %s'
  1357. % (M.shape, A.shape))
  1358. if np.dtype(M.dtype).char.lower() != np.dtype(A.dtype).char.lower():
  1359. warnings.warn('M does not have the same type precision as A. '
  1360. 'This may adversely affect ARPACK convergence')
  1361. n = A.shape[0]
  1362. if k <= 0:
  1363. raise ValueError("k must be greater than 0.")
  1364. if k >= n:
  1365. warnings.warn("k >= N for N * N square matrix. "
  1366. "Attempting to use scipy.linalg.eigh instead.",
  1367. RuntimeWarning)
  1368. if issparse(A):
  1369. raise TypeError("Cannot use scipy.linalg.eigh for sparse A with "
  1370. "k >= N. Use scipy.linalg.eigh(A.toarray()) or"
  1371. " reduce k.")
  1372. if isinstance(A, LinearOperator):
  1373. raise TypeError("Cannot use scipy.linalg.eigh for LinearOperator "
  1374. "A with k >= N.")
  1375. if isinstance(M, LinearOperator):
  1376. raise TypeError("Cannot use scipy.linalg.eigh for LinearOperator "
  1377. "M with k >= N.")
  1378. return eigh(A, b=M, eigvals_only=not return_eigenvectors)
  1379. if sigma is None:
  1380. A = _aslinearoperator_with_dtype(A)
  1381. matvec = A.matvec
  1382. if OPinv is not None:
  1383. raise ValueError("OPinv should not be specified "
  1384. "with sigma = None.")
  1385. if M is None:
  1386. #standard eigenvalue problem
  1387. mode = 1
  1388. M_matvec = None
  1389. Minv_matvec = None
  1390. if Minv is not None:
  1391. raise ValueError("Minv should not be "
  1392. "specified with M = None.")
  1393. else:
  1394. #general eigenvalue problem
  1395. mode = 2
  1396. if Minv is None:
  1397. Minv_matvec = get_inv_matvec(M, hermitian=True, tol=tol)
  1398. else:
  1399. Minv = _aslinearoperator_with_dtype(Minv)
  1400. Minv_matvec = Minv.matvec
  1401. M_matvec = _aslinearoperator_with_dtype(M).matvec
  1402. else:
  1403. # sigma is not None: shift-invert mode
  1404. if Minv is not None:
  1405. raise ValueError("Minv should not be specified when sigma is")
  1406. # normal mode
  1407. if mode == 'normal':
  1408. mode = 3
  1409. matvec = None
  1410. if OPinv is None:
  1411. Minv_matvec = get_OPinv_matvec(A, M, sigma,
  1412. hermitian=True, tol=tol)
  1413. else:
  1414. OPinv = _aslinearoperator_with_dtype(OPinv)
  1415. Minv_matvec = OPinv.matvec
  1416. if M is None:
  1417. M_matvec = None
  1418. else:
  1419. M = _aslinearoperator_with_dtype(M)
  1420. M_matvec = M.matvec
  1421. # buckling mode
  1422. elif mode == 'buckling':
  1423. mode = 4
  1424. if OPinv is None:
  1425. Minv_matvec = get_OPinv_matvec(A, M, sigma,
  1426. hermitian=True, tol=tol)
  1427. else:
  1428. Minv_matvec = _aslinearoperator_with_dtype(OPinv).matvec
  1429. matvec = _aslinearoperator_with_dtype(A).matvec
  1430. M_matvec = None
  1431. # cayley-transform mode
  1432. elif mode == 'cayley':
  1433. mode = 5
  1434. matvec = _aslinearoperator_with_dtype(A).matvec
  1435. if OPinv is None:
  1436. Minv_matvec = get_OPinv_matvec(A, M, sigma,
  1437. hermitian=True, tol=tol)
  1438. else:
  1439. Minv_matvec = _aslinearoperator_with_dtype(OPinv).matvec
  1440. if M is None:
  1441. M_matvec = None
  1442. else:
  1443. M_matvec = _aslinearoperator_with_dtype(M).matvec
  1444. # unrecognized mode
  1445. else:
  1446. raise ValueError("unrecognized mode '%s'" % mode)
  1447. params = _SymmetricArpackParams(n, k, A.dtype.char, matvec, mode,
  1448. M_matvec, Minv_matvec, sigma,
  1449. ncv, v0, maxiter, which, tol)
  1450. with _ARPACK_LOCK:
  1451. while not params.converged:
  1452. params.iterate()
  1453. return params.extract(return_eigenvectors)
  1454. def _augmented_orthonormal_cols(x, k):
  1455. # extract the shape of the x array
  1456. n, m = x.shape
  1457. # create the expanded array and copy x into it
  1458. y = np.empty((n, m+k), dtype=x.dtype)
  1459. y[:, :m] = x
  1460. # do some modified gram schmidt to add k random orthonormal vectors
  1461. for i in range(k):
  1462. # sample a random initial vector
  1463. v = np.random.randn(n)
  1464. if np.iscomplexobj(x):
  1465. v = v + 1j*np.random.randn(n)
  1466. # subtract projections onto the existing unit length vectors
  1467. for j in range(m+i):
  1468. u = y[:, j]
  1469. v -= (np.dot(v, u.conj()) / np.dot(u, u.conj())) * u
  1470. # normalize v
  1471. v /= np.sqrt(np.dot(v, v.conj()))
  1472. # add v into the output array
  1473. y[:, m+i] = v
  1474. # return the expanded array
  1475. return y
  1476. def _augmented_orthonormal_rows(x, k):
  1477. return _augmented_orthonormal_cols(x.T, k).T
  1478. def _herm(x):
  1479. return x.T.conj()
  1480. def svds(A, k=6, ncv=None, tol=0, which='LM', v0=None,
  1481. maxiter=None, return_singular_vectors=True):
  1482. """Compute the largest k singular values/vectors for a sparse matrix.
  1483. Parameters
  1484. ----------
  1485. A : {sparse matrix, LinearOperator}
  1486. Array to compute the SVD on, of shape (M, N)
  1487. k : int, optional
  1488. Number of singular values and vectors to compute.
  1489. Must be 1 <= k < min(A.shape).
  1490. ncv : int, optional
  1491. The number of Lanczos vectors generated
  1492. ncv must be greater than k+1 and smaller than n;
  1493. it is recommended that ncv > 2*k
  1494. Default: ``min(n, max(2*k + 1, 20))``
  1495. tol : float, optional
  1496. Tolerance for singular values. Zero (default) means machine precision.
  1497. which : str, ['LM' | 'SM'], optional
  1498. Which `k` singular values to find:
  1499. - 'LM' : largest singular values
  1500. - 'SM' : smallest singular values
  1501. .. versionadded:: 0.12.0
  1502. v0 : ndarray, optional
  1503. Starting vector for iteration, of length min(A.shape). Should be an
  1504. (approximate) left singular vector if N > M and a right singular
  1505. vector otherwise.
  1506. Default: random
  1507. .. versionadded:: 0.12.0
  1508. maxiter : int, optional
  1509. Maximum number of iterations.
  1510. .. versionadded:: 0.12.0
  1511. return_singular_vectors : bool or str, optional
  1512. - True: return singular vectors (True) in addition to singular values.
  1513. .. versionadded:: 0.12.0
  1514. - "u": only return the u matrix, without computing vh (if N > M).
  1515. - "vh": only return the vh matrix, without computing u (if N <= M).
  1516. .. versionadded:: 0.16.0
  1517. Returns
  1518. -------
  1519. u : ndarray, shape=(M, k)
  1520. Unitary matrix having left singular vectors as columns.
  1521. If `return_singular_vectors` is "vh", this variable is not computed,
  1522. and None is returned instead.
  1523. s : ndarray, shape=(k,)
  1524. The singular values.
  1525. vt : ndarray, shape=(k, N)
  1526. Unitary matrix having right singular vectors as rows.
  1527. If `return_singular_vectors` is "u", this variable is not computed,
  1528. and None is returned instead.
  1529. Notes
  1530. -----
  1531. This is a naive implementation using ARPACK as an eigensolver
  1532. on A.H * A or A * A.H, depending on which one is more efficient.
  1533. Examples
  1534. --------
  1535. >>> from scipy.sparse import csc_matrix
  1536. >>> from scipy.sparse.linalg import svds, eigs
  1537. >>> A = csc_matrix([[1, 0, 0], [5, 0, 2], [0, -1, 0], [0, 0, 3]], dtype=float)
  1538. >>> u, s, vt = svds(A, k=2)
  1539. >>> s
  1540. array([ 2.75193379, 5.6059665 ])
  1541. >>> np.sqrt(eigs(A.dot(A.T), k=2)[0]).real
  1542. array([ 5.6059665 , 2.75193379])
  1543. """
  1544. if not (isinstance(A, LinearOperator) or isspmatrix(A)):
  1545. A = np.asarray(A)
  1546. n, m = A.shape
  1547. if k <= 0 or k >= min(n, m):
  1548. raise ValueError("k must be between 1 and min(A.shape), k=%d" % k)
  1549. if isinstance(A, LinearOperator):
  1550. if n > m:
  1551. X_dot = A.matvec
  1552. X_matmat = A.matmat
  1553. XH_dot = A.rmatvec
  1554. else:
  1555. X_dot = A.rmatvec
  1556. XH_dot = A.matvec
  1557. dtype = getattr(A, 'dtype', None)
  1558. if dtype is None:
  1559. dtype = A.dot(np.zeros([m,1])).dtype
  1560. # A^H * V; works around lack of LinearOperator.adjoint.
  1561. # XXX This can be slow!
  1562. def X_matmat(V):
  1563. out = np.empty((V.shape[1], m), dtype=dtype)
  1564. for i, col in enumerate(V.T):
  1565. out[i, :] = A.rmatvec(col.reshape(-1, 1)).T
  1566. return out.T
  1567. else:
  1568. if n > m:
  1569. X_dot = X_matmat = A.dot
  1570. XH_dot = _herm(A).dot
  1571. else:
  1572. XH_dot = A.dot
  1573. X_dot = X_matmat = _herm(A).dot
  1574. def matvec_XH_X(x):
  1575. return XH_dot(X_dot(x))
  1576. XH_X = LinearOperator(matvec=matvec_XH_X, dtype=A.dtype,
  1577. shape=(min(A.shape), min(A.shape)))
  1578. # Get a low rank approximation of the implicitly defined gramian matrix.
  1579. # This is not a stable way to approach the problem.
  1580. eigvals, eigvec = eigsh(XH_X, k=k, tol=tol ** 2, maxiter=maxiter,
  1581. ncv=ncv, which=which, v0=v0)
  1582. # In 'LM' mode try to be clever about small eigenvalues.
  1583. # Otherwise in 'SM' mode do not try to be clever.
  1584. if which == 'LM':
  1585. # Gramian matrices have real non-negative eigenvalues.
  1586. eigvals = np.maximum(eigvals.real, 0)
  1587. # Use the sophisticated detection of small eigenvalues from pinvh.
  1588. t = eigvec.dtype.char.lower()
  1589. factor = {'f': 1E3, 'd': 1E6}
  1590. cond = factor[t] * np.finfo(t).eps
  1591. cutoff = cond * np.max(eigvals)
  1592. # Get a mask indicating which eigenpairs are not degenerately tiny,
  1593. # and create the re-ordered array of thresholded singular values.
  1594. above_cutoff = (eigvals > cutoff)
  1595. nlarge = above_cutoff.sum()
  1596. nsmall = k - nlarge
  1597. slarge = np.sqrt(eigvals[above_cutoff])
  1598. s = np.zeros_like(eigvals)
  1599. s[:nlarge] = slarge
  1600. if not return_singular_vectors:
  1601. return s
  1602. if n > m:
  1603. vlarge = eigvec[:, above_cutoff]
  1604. ularge = X_matmat(vlarge) / slarge if return_singular_vectors != 'vh' else None
  1605. vhlarge = _herm(vlarge)
  1606. else:
  1607. ularge = eigvec[:, above_cutoff]
  1608. vhlarge = _herm(X_matmat(ularge) / slarge) if return_singular_vectors != 'u' else None
  1609. u = _augmented_orthonormal_cols(ularge, nsmall) if ularge is not None else None
  1610. vh = _augmented_orthonormal_rows(vhlarge, nsmall) if vhlarge is not None else None
  1611. elif which == 'SM':
  1612. s = np.sqrt(eigvals)
  1613. if not return_singular_vectors:
  1614. return s
  1615. if n > m:
  1616. v = eigvec
  1617. u = X_matmat(v) / s if return_singular_vectors != 'vh' else None
  1618. vh = _herm(v)
  1619. else:
  1620. u = eigvec
  1621. vh = _herm(X_matmat(u) / s) if return_singular_vectors != 'u' else None
  1622. else:
  1623. raise ValueError("which must be either 'LM' or 'SM'.")
  1624. return u, s, vh