lanczos
lanczos(
A,
v0=None,
deg=None,
rtol=1e-08,
orth=0,
sparse_mat=False,
return_basis=False,
seed=None,
dtype=None,
**kwargs,
)Lanczos method for symmetric tridiagonalization.
This function implements Paiges A27 variant (1) of the Lanczos method for tridiagonalizing linear operators, which builds a tridiagonal T from a symmetric A via an orthogonal change-of-basis Q:
\begin{align*}
K &= [v, Av, A^2 v, ..., A^{n-1}v] && \\
Q &= [q_1, q_2, ..., q_{n} ] \gets \mathrm{qr}(K) && \\
T &= Q^T A Q &&
\end{align*}
Note that unlike other Lanczos implementations (e.g. SciPy’s eigsh), which includes e.g. restarting, & deflation steps, this method simply executes deg steps of the Lanczos method with orth vector reorthogonalizations (per step) and returns the entries of T.
This implementation supports varying degrees of re-orthogonalization. In particular, orth=0 corresponds to no re-orthogonalization, orth < deg corresponds to partial re-orthogonalization, and orth >= deg corresponds to full re-orthogonalization. The number of matvecs scales linearly with deg and the number of inner-products scales quadratically with orth.
Parameters
| Name | Type | Description | Default |
|---|---|---|---|
| A | Union[np.ndarray, sparray, LinearOperator] | Symmetric operator to tridiagonalize. | required |
| v0 | Optional[np.ndarray] | Initial vector to orthogonalize against. | None |
| deg | Optional[int] | Size of the Krylov subspace to expand. | None |
| rtol | float | Relative tolerance to consider the invariant subspace as converged. | 1e-08 |
| orth | int | Number of additional Lanczos vectors to orthogonalize against. | 0 |
| sparse_mat | bool | Whether to output the tridiagonal matrix as a sparse matrix. | False |
| return_basis | bool | If True, returns the Krylov basis vectors Q. |
False |
| dtype | Optional[np.dtype] | The precision dtype to specialize the computation. | None |
Returns
| Name | Type | Description |
|---|---|---|
| tuple | A tuple (a,b) parameterizing the diagonal and off-diagonal of the tridiagonal Jacobi matrix. If return_basis=True, |
|
| tuple | the tuple (a,b), Q is returned, where Q represents an orthogonal basis for the degree-deg Krylov subspace. |
See Also
- scipy.linalg.eigh_tridiagonal : Eigenvalue solver for real symmetric tridiagonal matrices.
- operator.matrix_function : Approximates the action of a matrix function via the Lanczos method.
References
- Paige, Christopher C. “Computational variants of the Lanczos method for the eigenproblem.” IMA Journal of Applied Mathematics 10.3 (1972): 373-381.