bundles / scipy 1.17.1 / scipy / sparse / linalg / _eigen / lobpcg / lobpcg / lobpcg
function
scipy.sparse.linalg._eigen.lobpcg.lobpcg:lobpcg
Signature
def lobpcg ( A , X , B = None , M = None , Y = None , tol = None , maxiter = None , largest = True , verbosityLevel = 0 , retLambdaHistory = False , retResidualNormsHistory = False , restartControl = 20 ) Summary
Locally Optimal Block Preconditioned Conjugate Gradient Method (LOBPCG).
Extended Summary
LOBPCG is a preconditioned eigensolver for large real symmetric and complex Hermitian definite generalized eigenproblems.
Parameters
A: {sparse matrix, ndarray, LinearOperator, callable object}The Hermitian linear operator of the problem, usually given by a sparse matrix. Often called the "stiffness matrix".
X: ndarray, float32 or float64Initial approximation to the
keigenvectors (non-sparse). IfAhasshape=(n,n)thenXmust haveshape=(n,k).B: {sparse matrix, ndarray, LinearOperator, callable object}Optional. By default
B = None, which is equivalent to identity. The right hand side operator in a generalized eigenproblem if present. Often called the "mass matrix". Must be Hermitian positive definite.M: {sparse matrix, ndarray, LinearOperator, callable object}Optional. By default
M = None, which is equivalent to identity. Preconditioner aiming to accelerate convergence.Y: ndarray, float32 or float64, default: NoneAn
n-by-sizeYndarray of constraints withsizeY < n. The iterations will be performed in theB-orthogonal complement of the column-space ofY.Ymust be full rank if present.tol: scalar, optionalThe default is
tol=n*sqrt(eps). Solver tolerance for the stopping criterion.maxiter: int, default: 20Maximum number of iterations.
largest: bool, default: TrueWhen True, solve for the largest eigenvalues, otherwise the smallest.
verbosityLevel: int, optionalBy default
verbosityLevel=0no output. Controls the solver standard/screen output.retLambdaHistory: bool, default: FalseWhether to return iterative eigenvalue history.
retResidualNormsHistory: bool, default: FalseWhether to return iterative history of residual norms.
restartControl: int, optional.Iterations restart if the residuals jump
2**restartControltimes compared to the smallest recorded inretResidualNormsHistory. The default isrestartControl=20, making the restarts rare for backward compatibility.
Returns
lambda: ndarray of the shape ``(k, )``.Array of
kapproximate eigenvalues.v: ndarray of the same shape as ``X.shape``.An array of
kapproximate eigenvectors.lambdaHistory: ndarray, optional.The eigenvalue history, if
retLambdaHistoryisTrue.ResidualNormsHistory: ndarray, optional.The history of residual norms, if
retResidualNormsHistoryisTrue.
Notes
The iterative loop runs maxit=maxiter (20 if maxit=None) iterations at most and finishes earlier if the tolerance is met. Breaking backward compatibility with the previous version, LOBPCG now returns the block of iterative vectors with the best accuracy rather than the last one iterated, as a cure for possible divergence.
If X.dtype == np.float32 and user-provided operations/multiplications by A, B, and M all preserve the np.float32 data type, all the calculations and the output are in np.float32.
The size of the iteration history output equals to the number of the best (limited by maxit) iterations plus 3: initial, final, and postprocessing.
If both retLambdaHistory and retResidualNormsHistory are True, the return tuple has the following format (lambda, V, lambda history, residual norms history).
In the following n denotes the matrix size and k the number of required eigenvalues (smallest or largest).
The LOBPCG code internally solves eigenproblems of the size 3k on every iteration by calling the dense eigensolver eigh, so if k is not small enough compared to n, it makes no sense to call the LOBPCG code. Moreover, if one calls the LOBPCG algorithm for 5k > n, it would likely break internally, so the code calls the standard function eigh instead. It is not that n should be large for the LOBPCG to work, but rather the ratio n / k should be large. It you call LOBPCG with k=1 and n=10, it works though n is small. The method is intended for extremely large n / k.
The convergence speed depends basically on three factors:
Quality of the initial approximations
Xto the seeking eigenvectors. Randomly distributed around the origin vectors work well if no better choice is known.Relative separation of the desired eigenvalues from the rest of the eigenvalues. One can vary
kto improve the separation.Proper preconditioning to shrink the spectral spread. For example, a rod vibration test problem (under tests directory) is ill-conditioned for large
n, so convergence will be slow, unless efficient preconditioning is used. For this specific problem, a good simple preconditioner function would be a linear solve forA, which is easy to code sinceAis tridiagonal.
Examples
Our first example is minimalistic - find the largest eigenvalue of a diagonal matrix by solving the non-generalized eigenvalue problem ``A x = lambda x`` without constraints or preconditioning.import numpy as np from scipy.sparse import diags_array from scipy.sparse.linalg import LinearOperator, aslinearoperator from scipy.sparse.linalg import lobpcg✓
n = 100
✓vals = np.arange(1, n + 1).astype(np.int16)
✓A = diags_array(vals, offsets=0, shape=(n, n), dtype=None) A.toarray()✓
k = 1 rng = np.random.default_rng() X = rng.normal(size=(n, k)) X = X.astype(np.float32)✓
eigenvalues, _ = lobpcg(A, X, maxiter=60)
✓eigenvalues
✗A_lambda = lambda X: vals[:, np.newaxis] * X
✓def A_matmat(X): return vals[:, np.newaxis] * X✓
eigenvalues, _ = lobpcg(A_lambda, X, maxiter=60)
✓eigenvalues
✗eigenvalues, _ = lobpcg(A_matmat, X, maxiter=60)
✓eigenvalues
✗A_lo = LinearOperator((n, n), matvec=A_matmat, matmat=A_matmat, dtype=np.int16) eigenvalues, _ = lobpcg(A_lo, X, maxiter=80)✓
eigenvalues
✗eigenvalues, _ = lobpcg(aslinearoperator(A), X, maxiter=80)
✓eigenvalues
✗k = 3 X = np.random.default_rng().normal(size=(n, k))✓
eigenvalues, _ = lobpcg(A, X, largest=False, maxiter=90) print(eigenvalues)✓
Y = np.eye(n, 3)
✓inv_vals = 1./vals inv_vals = inv_vals.astype(np.float32) M = lambda X: inv_vals[:, np.newaxis] * X✓
eigenvalues, _ = lobpcg(A_matmat, X, Y=Y, largest=False, maxiter=80) eigenvalues eigenvalues.dtype✓
eigenvalues, _ = lobpcg(A_matmat, X, Y=Y, M=M, largest=False, maxiter=20) eigenvalues✓
vals = vals - 50 X = rng.normal(size=(n, k)) eigenvalues, _ = lobpcg(A_matmat, X, largest=False, maxiter=99)✓
eigenvalues
✗eigenvalues, _ = lobpcg(A_matmat, X, largest=True, maxiter=99)
✓eigenvalues
✗Aliases
-
scipy.sparse.linalg.lobpcg