low_rank_toolbox.krylov.spaces.extended_krylov_space
Extended Krylov subspace construction.
Author: Benjamin Carrel, University of Geneva, 2022-2023
Classes
|
Extended Krylov space. |
- class low_rank_toolbox.krylov.spaces.extended_krylov_space.ExtendedKrylovSpace(A, X, invA=None, **extra_args)[source]
Bases:
SpaceStructureExtended Krylov space.
Constructs the extended Krylov space by combining a standard Krylov space and an inverted Krylov space:
EK_m(A, X) = span{X, AX, A^2X, …, A^(m-1)X, A^(-1)X, A^(-2)X, …, A^(-m)X}
This space is particularly useful for problems requiring information from both the range of A and its inverse.
Algorithm Selection: When is_symmetric=True, both component Krylov spaces (standard and inverted) use the Lanczos algorithm, providing approximately 2x speedup over the non-symmetric case. The combined basis maintains orthogonality across both components.
How to use
Create an instance: EK = ExtendedKrylovSpace(A, X, invA=None)
Augment the basis as needed: EK.augment_basis()
Access the basis via EK.basis or EK.Q
- A
Sparse matrix of shape (n, n)
- Type:
spmatrix
- X
Initial vector or matrix of shape (n, r)
- Type:
ndarray
- invA
Function that computes A^(-1) * v for a given vector/matrix v
- Type:
callable
- krylov_space
The standard Krylov space component
- Type:
- inverted_krylov_space
The inverted Krylov space component
- Type:
- Q
Combined orthonormal basis of the extended Krylov space
- Type:
ndarray
- Q1
Basis of the standard Krylov space
- Type:
ndarray
- Q2
Basis of the inverted Krylov space
- Type:
ndarray
- basis
Pointer to Q
- Type:
ndarray
Examples
>>> import numpy as np >>> from scipy.sparse import csr_matrix >>> from low_rank_toolbox.krylov import ExtendedKrylovSpace >>> # Create a sparse matrix >>> A = csr_matrix([[4, 1], [1, 3]]) >>> X = np.array([[1.0], [0.0]]) >>> # Extended Krylov includes both A^k*X and A^(-k)*X >>> EK = ExtendedKrylovSpace(A, X) >>> EK.augment_basis() # Adds both A*X and A^(-1)*X >>> EK.size # Total size is sum of both components 4 >>> # Access individual components >>> EK.Q1.shape # Standard Krylov part (2, 2) >>> EK.Q2.shape # Inverted Krylov part (2, 2) >>> # Combined basis is orthonormal >>> np.allclose(EK.Q.T @ EK.Q, np.eye(EK.Q.shape[1])) True
- property H1: ndarray
Upper Hessenberg matrix from Arnoldi for the Krylov space.
- Returns:
The Hessenberg matrix (non-symmetric case only).
- Return type:
ndarray
- property H2: ndarray
Upper Hessenberg matrix from Arnoldi for the inverted Krylov space.
- Returns:
The Hessenberg matrix (non-symmetric case only).
- Return type:
ndarray
- property Q: ndarray
Combined orthonormal basis of the extended Krylov space.
Concatenates Q1 and Q2, then orthonormalizes the result. The result is cached and only recomputed when the space size changes.
- Returns:
Matrix of shape (n, size) containing the orthonormal basis.
- Return type:
ndarray
- property Q1: ndarray
Basis of the standard Krylov space component.
- Returns:
Matrix of shape (n, m*r) containing the Krylov basis vectors.
- Return type:
ndarray
- property Q2: ndarray
Basis of the inverted Krylov space component.
- Returns:
Matrix of shape (n, m*r) containing the inverted Krylov basis vectors.
- Return type:
ndarray
- __init__(A, X, invA=None, **extra_args)[source]
- Parameters:
A (spmatrix) – The matrix A of the linear system.
X (ndarray) – The basis of the Krylov space.
invA (callable) – The function that computes the action of A^(-1) on a vector, or a matrix.
- augment_basis()[source]
Augment the basis of the extended Krylov space.
Augments both the standard Krylov space and the inverted Krylov space by adding the next block of basis vectors to each component.
- Parameters:
A (spmatrix)
X (ndarray)
invA (Optional[Callable])