low_rank_toolbox.krylov.spaces
Krylov Subspace Constructions.
This submodule provides classes for constructing and working with various types of Krylov subspaces, which are fundamental for iterative matrix computations.
Classes
KrylovSpace : Standard Krylov subspace K_m(A, b) InvertedKrylovSpace : Inverted Krylov subspace K_m(A^{-1}, b) ExtendedKrylovSpace : Extended Krylov subspace EK_m(A, b) RationalKrylovSpace : Rational Krylov subspace with arbitrary poles
- class low_rank_toolbox.krylov.spaces.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])
- class low_rank_toolbox.krylov.spaces.InvertedKrylovSpace(A, X, invA=None, **extra_args)[source]
Bases:
KrylovSpaceInverted Krylov space.
- Constructs the inverted Krylov space:
IK_m(A, X) = span{A^(-1)X, A^(-2)X, …, A^(-m)X}
This class wraps the KrylovSpace class by providing a custom matvec function that computes A^(-1) * v instead of A * v.
How to use
Initialize the inverted Krylov space with matrix A and vector/matrix X.
Augment the basis as needed with the augment_basis method.
Access the basis via the basis or Q attribute.
- 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
- Q
Orthonormal 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 InvertedKrylovSpace >>> # Create a sparse matrix >>> A = csr_matrix([[4, 1, 0], [1, 3, 1], [0, 1, 2]]) >>> X = np.array([[1.0], [0.0], [0.0]]) >>> # Build inverted Krylov space >>> IK = InvertedKrylovSpace(A, X) >>> # Start with A^(-1)*X >>> IK.Q.shape (3, 1) >>> # Augment with A^(-2)*X >>> IK.augment_basis() >>> IK.Q.shape (3, 2) >>> # Basis is orthonormal >>> np.allclose(IK.Q.T @ IK.Q, np.eye(2)) True
- class low_rank_toolbox.krylov.spaces.KrylovSpace(A, X, **extra_args)[source]
Bases:
SpaceStructureClass for (block) Krylov spaces. The definition of a Krylov space of size $m$ is the following:
K_m(A,x) = span{x, A x, A^2 x, …, A^(m-1) x}
where $A$ is a sparse matrix, and $x$ is a vector or a matrix.
Algorithm Selection: - If A is symmetric (is_symmetric=True): Uses Lanczos algorithm with 3-term recurrence
Stores tridiagonal coefficients (alpha, beta) in O(m) space
Orthogonalization cost: O(n*r) per iteration
If A is non-symmetric (is_symmetric=False): Uses Arnoldi algorithm with full orthogonalization - Stores upper Hessenberg matrix H in O(m²) space - Orthogonalization cost: O(n*m*r) per iteration
Performance: For symmetric problems, Lanczos is faster and uses significantly less memory than Arnoldi.
How to use
Initialize the Krylov space with the matrix $A$ and the vector $x$.
Augment the basis of the Krylov space as needed with the method augment_basis.
The basis is stored in the attribute ‘basis’, or ‘Q’ for short.
- A
Matrix of shape (n,n)
- Type:
spmatrix
- X
Vector of shape (n,1) or (n,r)
- Type:
ndarray
- Q
Matrix of shape (n,m) or (n,m*r) containing the basis of the Krylov space
- Type:
ndarray
- basis
Pointer to Q
- Type:
ndarray
- _alpha
Lanczos diagonal coefficients
- Type:
ndarray (symmetric only)
- _beta
Lanczos off-diagonal coefficients
- Type:
ndarray (symmetric only)
- H
Upper Hessenberg matrix from Arnoldi
- Type:
ndarray (non-symmetric only)
Examples
>>> import numpy as np >>> from scipy.sparse import csr_matrix >>> from low_rank_toolbox.krylov import KrylovSpace >>> # Non-symmetric case (uses Arnoldi) >>> A = csr_matrix([[1, 2, 0], [0, 3, 1], [1, 0, 2]]) >>> x = np.array([1.0, 0.0, 0.0]) >>> K = KrylovSpace(A, x, is_symmetric=False) >>> K.augment_basis() # Add next basis vector >>> K.Q.shape (3, 2) >>> # Symmetric case (uses Lanczos - more efficient) >>> A_sym = csr_matrix([[4, 1, 0], [1, 3, 1], [0, 1, 2]]) >>> K_sym = KrylovSpace(A_sym, x, is_symmetric=True) >>> K_sym.augment_basis() >>> K_sym.Q.shape (3, 2)
- __init__(A, X, **extra_args)[source]
Initialize a Krylov Space where X is a vector or a matrix
- Parameters:
A (spmatrix) – Sparse matrix of shape (n,n)
X (ndarray) – Vector or matrix of shape (n,) or (n,r)
extra_args (dict) – Extra arguments
arguments (Extra)
---------------
is_symmetric (bool) – True if A is symmetric (uses Lanczos), False otherwise (uses Arnoldi). Lanczos is much cheaper for symmetric problems.
matvec (callable) – Function for the matrix-vector product
- Return type:
None
- augment_basis()[source]
Augment the basis of the Krylov space.
Adds the next block of r basis vectors to the Krylov space using: - Lanczos algorithm if A is symmetric - Arnoldi algorithm if A is non-symmetric
The new basis vectors are computed as A * Q[:, (m-1)*r:m*r] and then orthogonalized against the existing basis.
- Return type:
Notes
If the next basis would exceed the dimension of the matrix, a warning is issued and the method returns without modifying the basis.
- property basis
The orthonormal basis of the Krylov space.
- Returns:
Matrix of shape (n, m*r) containing the basis vectors.
- Return type:
ndarray
- Parameters:
A (spmatrix)
X (ndarray)
- class low_rank_toolbox.krylov.spaces.RationalKrylovSpace(A, X, poles, inverse_only=False, **extra_args)[source]
Bases:
SpaceStructureRational Krylov space.
- Constructs the rational Krylov space:
RK_m(A, X) = span{X, (A - p_1*I)^{-1}X, …, prod_{i=1}^{m-1}(A - p_i*I)^{-1}X}
where p_1, …, p_{m-1} are specified poles (shifts). This space generalizes the standard Krylov space by using rational functions of A instead of polynomials.
How to use
Initialize with matrix A, vector/matrix X, and list of poles.
Augment the basis as needed with augment_basis.
Access the basis via basis or Q.
- A
Sparse matrix of shape (n, n)
- Type:
spmatrix
- X
Initial vector or matrix of shape (n, r)
- Type:
ndarray
- poles
Array of poles (shifts) for the rational Krylov space
- Type:
ndarray
- Q
Orthonormal basis of shape (n, m*r)
- Type:
ndarray
- H
Upper triangular matrix from QR factorization
- 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 RationalKrylovSpace >>> # Create a sparse matrix >>> A = csr_matrix([[4, 1, 0], [1, 3, 1], [0, 1, 2]]) >>> X = np.array([[1.0], [0.0], [0.0]]) >>> # Choose poles to emphasize specific spectral regions >>> # (e.g., near eigenvalues of interest) >>> poles = [1.0, 2.0, 3.0] >>> RK = RationalKrylovSpace(A, X, poles) >>> # Each augmentation uses the next pole >>> RK.augment_basis() # Uses pole 1.0 >>> RK.Q.shape (3, 2) >>> RK.augment_basis() # Uses pole 2.0 >>> RK.Q.shape (3, 3) >>> # Verify orthogonality >>> np.allclose(RK.Q.T @ RK.Q, np.eye(3)) True
- __init__(A, X, poles, inverse_only=False, **extra_args)[source]
Initialize a rational Krylov space
- Parameters:
A (spmatrix) – Sparse matrix of shape (n, n)
X (ndarray) – Vector or matrix of shape (n, r)
poles (list) – Poles (shifts) for the rational Krylov space
inverse_only (bool, optional) – If True, solves (A - p*I) * v_new = v_old. If False, solves (A - p*I) * v_new = A * v_old. Default is False. The False case is useful for approximating matrix functions, while True is for solving linear systems.
extra_args (dict) – Extra arguments
arguments (Extra)
---------------
symmetric (bool) – True if A is symmetric, False otherwise (not used for rational Krylov)
- Return type:
None
- augment_basis()[source]
Augment the basis of the rational Krylov space.
- Adds the next block of r basis vectors by solving:
(A - p_{m-1}*I) * v_new = v_old
where p_{m-1} is the next pole in the sequence. The new vectors are then orthogonalized against the existing basis using QR factorization.
- Raises:
ValueError – If the space would exceed the problem dimension or if there are not enough poles specified.
Modules
Extended Krylov subspace construction. |
|
Inverted Krylov subspace construction. |
|
Standard Krylov subspace construction and operations. |
|
Rational Krylov subspace construction with arbitrary poles. |
|
Base class for space structures. |