arpack.py 66 KB

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