low_rank_toolbox.krylov.spaces.rational_krylov_space
Rational Krylov subspace construction with arbitrary poles.
Author: Benjamin Carrel, University of Geneva, 2022-2023
Classes
|
Rational Krylov space. |
- class low_rank_toolbox.krylov.spaces.rational_krylov_space.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.