cupyx.scipy.sparse.linalg.lobpcg#
- cupyx.scipy.sparse.linalg.lobpcg(A, X, B=None, M=None, Y=None, tol=None, maxiter=None, largest=True, verbosityLevel=0, retLambdaHistory=False, retResidualNormsHistory=False)[source]#
Locally Optimal Block Preconditioned Conjugate Gradient Method (LOBPCG)
LOBPCG is a preconditioned eigensolver for large symmetric positive definite (SPD) generalized eigenproblems.
- Parameters:
A (array-like) – The symmetric linear operator of the problem, usually a sparse matrix. Can be of the following types - cupy.ndarray - cupyx.scipy.sparse.csr_matrix - cupy.scipy.sparse.linalg.LinearOperator
X (cupy.ndarray) – Initial approximation to the
k
eigenvectors (non-sparse). If A hasshape=(n,n)
then X should have shapeshape=(n,k)
.B (array-like) – The right hand side operator in a generalized eigenproblem. By default,
B = Identity
. Can be of following types: - cupy.ndarray - cupyx.scipy.sparse.csr_matrix - cupy.scipy.sparse.linalg.LinearOperatorM (array-like) – Preconditioner to A; by default
M = Identity
. M should approximate the inverse of A. Can be of the following types: - cupy.ndarray - cupyx.scipy.sparse.csr_matrix - cupy.scipy.sparse.linalg.LinearOperatorY (cupy.ndarray) – n-by-sizeY matrix of constraints (non-sparse), sizeY < n The iterations will be performed in the B-orthogonal complement of the column-space of Y. Y must be full rank.
tol (float) – Solver tolerance (stopping criterion). The default is
tol=n*sqrt(eps)
.maxiter (int) – Maximum number of iterations. The default is
maxiter = 20
.largest (bool) – When True, solve for the largest eigenvalues, otherwise the smallest.
verbosityLevel (int) – Controls solver output. The default is
verbosityLevel=0
.retLambdaHistory (bool) – Whether to return eigenvalue history. Default is False.
retResidualNormsHistory (bool) – Whether to return history of residual norms. Default is False.
- Returns:
w (cupy.ndarray): Array of
k
eigenvaluesv (cupy.ndarray) An array of
k
eigenvectors. v has the same shape as X.lambdas (list of cupy.ndarray): The eigenvalue history, if retLambdaHistory is True.
rnorms (list of cupy.ndarray): The history of residual norms, if retResidualNormsHistory is True.
- Return type:
See also
Note
If both
retLambdaHistory
andretResidualNormsHistory
are True the return tuple has the following format(lambda, V, lambda history, residual norms history)
.