lobpcg.py 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582
  1. """
  2. Pure SciPy implementation of Locally Optimal Block Preconditioned Conjugate
  3. Gradient Method (LOBPCG), see
  4. https://bitbucket.org/joseroman/blopex
  5. License: BSD
  6. Authors: Robert Cimrman, Andrew Knyazev
  7. Examples in tests directory contributed by Nils Wagner.
  8. """
  9. from __future__ import division, print_function, absolute_import
  10. import numpy as np
  11. from numpy.testing import assert_allclose
  12. from scipy._lib.six import xrange
  13. from scipy.linalg import inv, eigh, cho_factor, cho_solve, cholesky
  14. from scipy.sparse.linalg import aslinearoperator, LinearOperator
  15. __all__ = ['lobpcg']
  16. def save(ar, fileName):
  17. # Used only when verbosity level > 10.
  18. from numpy import savetxt
  19. savetxt(fileName, ar)
  20. def _report_nonhermitian(M, a, b, name):
  21. """
  22. Report if `M` is not a hermitian matrix given the tolerances `a`, `b`.
  23. """
  24. from scipy.linalg import norm
  25. md = M - M.T.conj()
  26. nmd = norm(md, 1)
  27. tol = np.spacing(max(10**a, (10**b)*norm(M, 1)))
  28. if nmd > tol:
  29. print('matrix %s is not sufficiently Hermitian for a=%d, b=%d:'
  30. % (name, a, b))
  31. print('condition: %.e < %e' % (nmd, tol))
  32. ##
  33. # 21.05.2007, c
  34. def as2d(ar):
  35. """
  36. If the input array is 2D return it, if it is 1D, append a dimension,
  37. making it a column vector.
  38. """
  39. if ar.ndim == 2:
  40. return ar
  41. else: # Assume 1!
  42. aux = np.array(ar, copy=False)
  43. aux.shape = (ar.shape[0], 1)
  44. return aux
  45. def _makeOperator(operatorInput, expectedShape):
  46. """Takes a dense numpy array or a sparse matrix or
  47. a function and makes an operator performing matrix * blockvector
  48. products.
  49. Examples
  50. --------
  51. >>> A = _makeOperator( arrayA, (n, n) )
  52. >>> vectorB = A( vectorX )
  53. """
  54. if operatorInput is None:
  55. def ident(x):
  56. return x
  57. operator = LinearOperator(expectedShape, ident, matmat=ident)
  58. else:
  59. operator = aslinearoperator(operatorInput)
  60. if operator.shape != expectedShape:
  61. raise ValueError('operator has invalid shape')
  62. return operator
  63. def _applyConstraints(blockVectorV, factYBY, blockVectorBY, blockVectorY):
  64. """Changes blockVectorV in place."""
  65. gramYBV = np.dot(blockVectorBY.T.conj(), blockVectorV)
  66. tmp = cho_solve(factYBY, gramYBV)
  67. blockVectorV -= np.dot(blockVectorY, tmp)
  68. def _b_orthonormalize(B, blockVectorV, blockVectorBV=None, retInvR=False):
  69. if blockVectorBV is None:
  70. if B is not None:
  71. blockVectorBV = B(blockVectorV)
  72. else:
  73. blockVectorBV = blockVectorV # Shared data!!!
  74. gramVBV = np.dot(blockVectorV.T.conj(), blockVectorBV)
  75. gramVBV = cholesky(gramVBV)
  76. gramVBV = inv(gramVBV, overwrite_a=True)
  77. # gramVBV is now R^{-1}.
  78. blockVectorV = np.dot(blockVectorV, gramVBV)
  79. if B is not None:
  80. blockVectorBV = np.dot(blockVectorBV, gramVBV)
  81. if retInvR:
  82. return blockVectorV, blockVectorBV, gramVBV
  83. else:
  84. return blockVectorV, blockVectorBV
  85. def _get_indx(_lambda, num, largest):
  86. """Get `num` indices into `_lambda` depending on `largest` option."""
  87. ii = np.argsort(_lambda)
  88. if largest:
  89. ii = ii[:-num-1:-1]
  90. else:
  91. ii = ii[:num]
  92. return ii
  93. def lobpcg(A, X,
  94. B=None, M=None, Y=None,
  95. tol=None, maxiter=20,
  96. largest=True, verbosityLevel=0,
  97. retLambdaHistory=False, retResidualNormsHistory=False):
  98. """Locally Optimal Block Preconditioned Conjugate Gradient Method (LOBPCG)
  99. LOBPCG is a preconditioned eigensolver for large symmetric positive
  100. definite (SPD) generalized eigenproblems.
  101. Parameters
  102. ----------
  103. A : {sparse matrix, dense matrix, LinearOperator}
  104. The symmetric linear operator of the problem, usually a
  105. sparse matrix. Often called the "stiffness matrix".
  106. X : array_like
  107. Initial approximation to the k eigenvectors. If A has
  108. shape=(n,n) then X should have shape shape=(n,k).
  109. B : {dense matrix, sparse matrix, LinearOperator}, optional
  110. the right hand side operator in a generalized eigenproblem.
  111. by default, B = Identity
  112. often called the "mass matrix"
  113. M : {dense matrix, sparse matrix, LinearOperator}, optional
  114. preconditioner to A; by default M = Identity
  115. M should approximate the inverse of A
  116. Y : array_like, optional
  117. n-by-sizeY matrix of constraints, sizeY < n
  118. The iterations will be performed in the B-orthogonal complement
  119. of the column-space of Y. Y must be full rank.
  120. Returns
  121. -------
  122. w : array
  123. Array of k eigenvalues
  124. v : array
  125. An array of k eigenvectors. V has the same shape as X.
  126. Other Parameters
  127. ----------------
  128. tol : scalar, optional
  129. Solver tolerance (stopping criterion)
  130. by default: tol=n*sqrt(eps)
  131. maxiter : integer, optional
  132. maximum number of iterations
  133. by default: maxiter=min(n,20)
  134. largest : bool, optional
  135. when True, solve for the largest eigenvalues, otherwise the smallest
  136. verbosityLevel : integer, optional
  137. controls solver output. default: verbosityLevel = 0.
  138. retLambdaHistory : boolean, optional
  139. whether to return eigenvalue history
  140. retResidualNormsHistory : boolean, optional
  141. whether to return history of residual norms
  142. Examples
  143. --------
  144. Solve A x = lambda B x with constraints and preconditioning.
  145. >>> from scipy.sparse import spdiags, issparse
  146. >>> from scipy.sparse.linalg import lobpcg, LinearOperator
  147. >>> n = 100
  148. >>> vals = [np.arange(n, dtype=np.float64) + 1]
  149. >>> A = spdiags(vals, 0, n, n)
  150. >>> A.toarray()
  151. array([[ 1., 0., 0., ..., 0., 0., 0.],
  152. [ 0., 2., 0., ..., 0., 0., 0.],
  153. [ 0., 0., 3., ..., 0., 0., 0.],
  154. ...,
  155. [ 0., 0., 0., ..., 98., 0., 0.],
  156. [ 0., 0., 0., ..., 0., 99., 0.],
  157. [ 0., 0., 0., ..., 0., 0., 100.]])
  158. Constraints.
  159. >>> Y = np.eye(n, 3)
  160. Initial guess for eigenvectors, should have linearly independent
  161. columns. Column dimension = number of requested eigenvalues.
  162. >>> X = np.random.rand(n, 3)
  163. Preconditioner -- inverse of A (as an abstract linear operator).
  164. >>> invA = spdiags([1./vals[0]], 0, n, n)
  165. >>> def precond( x ):
  166. ... return invA * x
  167. >>> M = LinearOperator(matvec=precond, shape=(n, n), dtype=float)
  168. Here, ``invA`` could of course have been used directly as a preconditioner.
  169. Let us then solve the problem:
  170. >>> eigs, vecs = lobpcg(A, X, Y=Y, M=M, tol=1e-4, maxiter=40, largest=False)
  171. >>> eigs
  172. array([ 4., 5., 6.])
  173. Note that the vectors passed in Y are the eigenvectors of the 3 smallest
  174. eigenvalues. The results returned are orthogonal to those.
  175. Notes
  176. -----
  177. If both retLambdaHistory and retResidualNormsHistory are True,
  178. the return tuple has the following format
  179. (lambda, V, lambda history, residual norms history).
  180. In the following ``n`` denotes the matrix size and ``m`` the number
  181. of required eigenvalues (smallest or largest).
  182. The LOBPCG code internally solves eigenproblems of the size 3``m`` on every
  183. iteration by calling the "standard" dense eigensolver, so if ``m`` is not
  184. small enough compared to ``n``, it does not make sense to call the LOBPCG
  185. code, but rather one should use the "standard" eigensolver,
  186. e.g. numpy or scipy function in this case.
  187. If one calls the LOBPCG algorithm for 5``m``>``n``,
  188. it will most likely break internally, so the code tries to call the standard
  189. function instead.
  190. It is not that n should be large for the LOBPCG to work, but rather the
  191. ratio ``n``/``m`` should be large. It you call the LOBPCG code with ``m``=1
  192. and ``n``=10, it should work, though ``n`` is small. The method is intended
  193. for extremely large ``n``/``m``, see e.g., reference [28] in
  194. https://arxiv.org/abs/0705.2626
  195. The convergence speed depends basically on two factors:
  196. 1. How well relatively separated the seeking eigenvalues are
  197. from the rest of the eigenvalues.
  198. One can try to vary ``m`` to make this better.
  199. 2. How well conditioned the problem is. This can be changed by using proper
  200. preconditioning. For example, a rod vibration test problem (under tests
  201. directory) is ill-conditioned for large ``n``, so convergence will be
  202. slow, unless efficient preconditioning is used.
  203. For this specific problem, a good simple preconditioner function would
  204. be a linear solve for A, which is easy to code since A is tridiagonal.
  205. *Acknowledgements*
  206. lobpcg.py code was written by Robert Cimrman.
  207. Many thanks belong to Andrew Knyazev, the author of the algorithm,
  208. for lots of advice and support.
  209. References
  210. ----------
  211. .. [1] A. V. Knyazev (2001),
  212. Toward the Optimal Preconditioned Eigensolver: Locally Optimal
  213. Block Preconditioned Conjugate Gradient Method.
  214. SIAM Journal on Scientific Computing 23, no. 2,
  215. pp. 517-541. :doi:`10.1137/S1064827500366124`
  216. .. [2] A. V. Knyazev, I. Lashuk, M. E. Argentati, and E. Ovchinnikov (2007),
  217. Block Locally Optimal Preconditioned Eigenvalue Xolvers (BLOPEX)
  218. in hypre and PETSc. https://arxiv.org/abs/0705.2626
  219. .. [3] A. V. Knyazev's C and MATLAB implementations:
  220. https://bitbucket.org/joseroman/blopex
  221. """
  222. blockVectorX = X
  223. blockVectorY = Y
  224. residualTolerance = tol
  225. maxIterations = maxiter
  226. if blockVectorY is not None:
  227. sizeY = blockVectorY.shape[1]
  228. else:
  229. sizeY = 0
  230. # Block size.
  231. if len(blockVectorX.shape) != 2:
  232. raise ValueError('expected rank-2 array for argument X')
  233. n, sizeX = blockVectorX.shape
  234. if sizeX > n:
  235. raise ValueError('X column dimension exceeds the row dimension')
  236. A = _makeOperator(A, (n,n))
  237. B = _makeOperator(B, (n,n))
  238. M = _makeOperator(M, (n,n))
  239. if (n - sizeY) < (5 * sizeX):
  240. # warn('The problem size is small compared to the block size.' \
  241. # ' Using dense eigensolver instead of LOBPCG.')
  242. if blockVectorY is not None:
  243. raise NotImplementedError('The dense eigensolver '
  244. 'does not support constraints.')
  245. # Define the closed range of indices of eigenvalues to return.
  246. if largest:
  247. eigvals = (n - sizeX, n-1)
  248. else:
  249. eigvals = (0, sizeX-1)
  250. A_dense = A(np.eye(n))
  251. B_dense = None if B is None else B(np.eye(n))
  252. vals, vecs = eigh(A_dense, B_dense, eigvals=eigvals, check_finite=False)
  253. if largest:
  254. # Reverse order to be compatible with eigs() in 'LM' mode.
  255. vals = vals[::-1]
  256. vecs = vecs[:, ::-1]
  257. return vals, vecs
  258. if residualTolerance is None:
  259. residualTolerance = np.sqrt(1e-15) * n
  260. maxIterations = min(n, maxIterations)
  261. if verbosityLevel:
  262. aux = "Solving "
  263. if B is None:
  264. aux += "standard"
  265. else:
  266. aux += "generalized"
  267. aux += " eigenvalue problem with"
  268. if M is None:
  269. aux += "out"
  270. aux += " preconditioning\n\n"
  271. aux += "matrix size %d\n" % n
  272. aux += "block size %d\n\n" % sizeX
  273. if blockVectorY is None:
  274. aux += "No constraints\n\n"
  275. else:
  276. if sizeY > 1:
  277. aux += "%d constraints\n\n" % sizeY
  278. else:
  279. aux += "%d constraint\n\n" % sizeY
  280. print(aux)
  281. ##
  282. # Apply constraints to X.
  283. if blockVectorY is not None:
  284. if B is not None:
  285. blockVectorBY = B(blockVectorY)
  286. else:
  287. blockVectorBY = blockVectorY
  288. # gramYBY is a dense array.
  289. gramYBY = np.dot(blockVectorY.T.conj(), blockVectorBY)
  290. try:
  291. # gramYBY is a Cholesky factor from now on...
  292. gramYBY = cho_factor(gramYBY)
  293. except Exception:
  294. raise ValueError('cannot handle linearly dependent constraints')
  295. _applyConstraints(blockVectorX, gramYBY, blockVectorBY, blockVectorY)
  296. ##
  297. # B-orthonormalize X.
  298. blockVectorX, blockVectorBX = _b_orthonormalize(B, blockVectorX)
  299. ##
  300. # Compute the initial Ritz vectors: solve the eigenproblem.
  301. blockVectorAX = A(blockVectorX)
  302. gramXAX = np.dot(blockVectorX.T.conj(), blockVectorAX)
  303. _lambda, eigBlockVector = eigh(gramXAX, check_finite=False)
  304. ii = _get_indx(_lambda, sizeX, largest)
  305. _lambda = _lambda[ii]
  306. eigBlockVector = np.asarray(eigBlockVector[:,ii])
  307. blockVectorX = np.dot(blockVectorX, eigBlockVector)
  308. blockVectorAX = np.dot(blockVectorAX, eigBlockVector)
  309. if B is not None:
  310. blockVectorBX = np.dot(blockVectorBX, eigBlockVector)
  311. ##
  312. # Active index set.
  313. activeMask = np.ones((sizeX,), dtype=bool)
  314. lambdaHistory = [_lambda]
  315. residualNormsHistory = []
  316. previousBlockSize = sizeX
  317. ident = np.eye(sizeX, dtype=A.dtype)
  318. ident0 = np.eye(sizeX, dtype=A.dtype)
  319. ##
  320. # Main iteration loop.
  321. blockVectorP = None # set during iteration
  322. blockVectorAP = None
  323. blockVectorBP = None
  324. for iterationNumber in xrange(maxIterations):
  325. if verbosityLevel > 0:
  326. print('iteration %d' % iterationNumber)
  327. aux = blockVectorBX * _lambda[np.newaxis,:]
  328. blockVectorR = blockVectorAX - aux
  329. aux = np.sum(blockVectorR.conjugate() * blockVectorR, 0)
  330. residualNorms = np.sqrt(aux)
  331. residualNormsHistory.append(residualNorms)
  332. ii = np.where(residualNorms > residualTolerance, True, False)
  333. activeMask = activeMask & ii
  334. if verbosityLevel > 2:
  335. print(activeMask)
  336. currentBlockSize = activeMask.sum()
  337. if currentBlockSize != previousBlockSize:
  338. previousBlockSize = currentBlockSize
  339. ident = np.eye(currentBlockSize, dtype=A.dtype)
  340. if currentBlockSize == 0:
  341. break
  342. if verbosityLevel > 0:
  343. print('current block size:', currentBlockSize)
  344. print('eigenvalue:', _lambda)
  345. print('residual norms:', residualNorms)
  346. if verbosityLevel > 10:
  347. print(eigBlockVector)
  348. activeBlockVectorR = as2d(blockVectorR[:,activeMask])
  349. if iterationNumber > 0:
  350. activeBlockVectorP = as2d(blockVectorP[:,activeMask])
  351. activeBlockVectorAP = as2d(blockVectorAP[:,activeMask])
  352. activeBlockVectorBP = as2d(blockVectorBP[:,activeMask])
  353. if M is not None:
  354. # Apply preconditioner T to the active residuals.
  355. activeBlockVectorR = M(activeBlockVectorR)
  356. ##
  357. # Apply constraints to the preconditioned residuals.
  358. if blockVectorY is not None:
  359. _applyConstraints(activeBlockVectorR,
  360. gramYBY, blockVectorBY, blockVectorY)
  361. ##
  362. # B-orthonormalize the preconditioned residuals.
  363. aux = _b_orthonormalize(B, activeBlockVectorR)
  364. activeBlockVectorR, activeBlockVectorBR = aux
  365. activeBlockVectorAR = A(activeBlockVectorR)
  366. if iterationNumber > 0:
  367. aux = _b_orthonormalize(B, activeBlockVectorP,
  368. activeBlockVectorBP, retInvR=True)
  369. activeBlockVectorP, activeBlockVectorBP, invR = aux
  370. activeBlockVectorAP = np.dot(activeBlockVectorAP, invR)
  371. ##
  372. # Perform the Rayleigh Ritz Procedure:
  373. # Compute symmetric Gram matrices:
  374. xaw = np.dot(blockVectorX.T.conj(), activeBlockVectorAR)
  375. waw = np.dot(activeBlockVectorR.T.conj(), activeBlockVectorAR)
  376. xbw = np.dot(blockVectorX.T.conj(), activeBlockVectorBR)
  377. if iterationNumber > 0:
  378. xap = np.dot(blockVectorX.T.conj(), activeBlockVectorAP)
  379. wap = np.dot(activeBlockVectorR.T.conj(), activeBlockVectorAP)
  380. pap = np.dot(activeBlockVectorP.T.conj(), activeBlockVectorAP)
  381. xbp = np.dot(blockVectorX.T.conj(), activeBlockVectorBP)
  382. wbp = np.dot(activeBlockVectorR.T.conj(), activeBlockVectorBP)
  383. gramA = np.bmat([[np.diag(_lambda), xaw, xap],
  384. [xaw.T.conj(), waw, wap],
  385. [xap.T.conj(), wap.T.conj(), pap]])
  386. gramB = np.bmat([[ident0, xbw, xbp],
  387. [xbw.T.conj(), ident, wbp],
  388. [xbp.T.conj(), wbp.T.conj(), ident]])
  389. else:
  390. gramA = np.bmat([[np.diag(_lambda), xaw],
  391. [xaw.T.conj(), waw]])
  392. gramB = np.bmat([[ident0, xbw],
  393. [xbw.T.conj(), ident]])
  394. if verbosityLevel > 0:
  395. _report_nonhermitian(gramA, 3, -1, 'gramA')
  396. _report_nonhermitian(gramB, 3, -1, 'gramB')
  397. if verbosityLevel > 10:
  398. save(gramA, 'gramA')
  399. save(gramB, 'gramB')
  400. # Solve the generalized eigenvalue problem.
  401. _lambda, eigBlockVector = eigh(gramA, gramB, check_finite=False)
  402. ii = _get_indx(_lambda, sizeX, largest)
  403. if verbosityLevel > 10:
  404. print(ii)
  405. _lambda = _lambda[ii]
  406. eigBlockVector = eigBlockVector[:,ii]
  407. lambdaHistory.append(_lambda)
  408. if verbosityLevel > 10:
  409. print('lambda:', _lambda)
  410. ## # Normalize eigenvectors!
  411. ## aux = np.sum( eigBlockVector.conjugate() * eigBlockVector, 0 )
  412. ## eigVecNorms = np.sqrt( aux )
  413. ## eigBlockVector = eigBlockVector / eigVecNorms[np.newaxis,:]
  414. # eigBlockVector, aux = _b_orthonormalize( B, eigBlockVector )
  415. if verbosityLevel > 10:
  416. print(eigBlockVector)
  417. ##
  418. # Compute Ritz vectors.
  419. if iterationNumber > 0:
  420. eigBlockVectorX = eigBlockVector[:sizeX]
  421. eigBlockVectorR = eigBlockVector[sizeX:sizeX+currentBlockSize]
  422. eigBlockVectorP = eigBlockVector[sizeX+currentBlockSize:]
  423. pp = np.dot(activeBlockVectorR, eigBlockVectorR)
  424. pp += np.dot(activeBlockVectorP, eigBlockVectorP)
  425. app = np.dot(activeBlockVectorAR, eigBlockVectorR)
  426. app += np.dot(activeBlockVectorAP, eigBlockVectorP)
  427. bpp = np.dot(activeBlockVectorBR, eigBlockVectorR)
  428. bpp += np.dot(activeBlockVectorBP, eigBlockVectorP)
  429. else:
  430. eigBlockVectorX = eigBlockVector[:sizeX]
  431. eigBlockVectorR = eigBlockVector[sizeX:]
  432. pp = np.dot(activeBlockVectorR, eigBlockVectorR)
  433. app = np.dot(activeBlockVectorAR, eigBlockVectorR)
  434. bpp = np.dot(activeBlockVectorBR, eigBlockVectorR)
  435. if verbosityLevel > 10:
  436. print(pp)
  437. print(app)
  438. print(bpp)
  439. blockVectorX = np.dot(blockVectorX, eigBlockVectorX) + pp
  440. blockVectorAX = np.dot(blockVectorAX, eigBlockVectorX) + app
  441. blockVectorBX = np.dot(blockVectorBX, eigBlockVectorX) + bpp
  442. blockVectorP, blockVectorAP, blockVectorBP = pp, app, bpp
  443. aux = blockVectorBX * _lambda[np.newaxis,:]
  444. blockVectorR = blockVectorAX - aux
  445. aux = np.sum(blockVectorR.conjugate() * blockVectorR, 0)
  446. residualNorms = np.sqrt(aux)
  447. if verbosityLevel > 0:
  448. print('final eigenvalue:', _lambda)
  449. print('final residual norms:', residualNorms)
  450. if retLambdaHistory:
  451. if retResidualNormsHistory:
  452. return _lambda, blockVectorX, lambdaHistory, residualNormsHistory
  453. else:
  454. return _lambda, blockVectorX, lambdaHistory
  455. else:
  456. if retResidualNormsHistory:
  457. return _lambda, blockVectorX, residualNormsHistory
  458. else:
  459. return _lambda, blockVectorX