pymor.algorithms.chol_qr
¶
Module Contents¶
- pymor.algorithms.chol_qr.shifted_chol_qr(A, product=None, return_R=False, maxiter=3, offset=0, orth_tol=None, recompute_shift=False, check_finite=True, copy=True, product_norm=None)[source]¶
Orthonormalize a
VectorArray
using the shifted CholeskyQR algorithm.This method computes a QR decomposition of a
VectorArray
via Cholesky factorizations of its Gramian matrix according to [FKN+20]. For ill-conditioned matrices, the Cholesky factorization will break down. In this case a diagonal shift will be applied to the Gramian.shifted_chol_qr(A, maxiter=3, orth_tol=None)
is equivalent to the shifted CholeskyQR3 algorithm (Algorithm 4.2 in [FKN+20]).shifted_chol_qr(A, maxiter=np.inf, orth_tol=<some_number>)
is equivalent to the shifted CholeskyQR algorithm (Algorithm 4.1 in [FKN+20]).shifted_chol_qr(A, product=<some_product>, maxiter=3, orth_tol=None)
is equivalent to the shifted CholeskyQR3 algorithm in an oblique inner product (Algorithm 5.1 in [FKN+20]).
Note
Note that the unshifted single iteration CholeskyQR algorithm is unstable. Setting
maxiter=1
will perform the simple shifted CholeskyQR algorithm (Algorithm 2.1 in [FKN+20]) which is stable but leads to poor orthogonality in general, i.e.\[\lVert Q^TQ - I \rVert_2 < 2,\](see Lemma 3.1 in [FKN+20]).
- Parameters:
A – The
VectorArray
which is to be orthonormalized.product – The inner product
Operator
w.r.t. which to orthonormalize. IfNone
, the Euclidean product is used.return_R – If
True
, the R matrix from QR decomposition is returned.maxiter – Maximum number of iterations. Defaults to 3.
offset – Assume that the first
offset
vectors are already orthonormal and apply the algorithm for the vectors starting atoffset + 1
.orth_tol – If not
None
, check if the resultingVectorArray
is really orthornormal and repeat the algorithm until the check passes ormaxiter
is reached.recompute_shift – If
False
, the shift is computed just once if the Cholesky decomposition fails and reused in possible further iterations. However, the initial shift might be too large for further iterations, which would lead to a non-orthonormal basis. IfTrue
, the shift is recomputed in iterations in which the Cholesky decomposition fails. Even for an ill-conditionedA
(at least for matrix condition numbers up to 10^20) is it able to compute an orthonormal basis at the cost of higher runtimes.check_finite – This argument is passed down to
SciPy linalg
functions. Disabling may give a performance gain, but may result in problems (crashes, non-termination) if the inputs do contain infinities or NaNs.copy – If
True
, create a copy ofA
instead of modifyingA
in-place.product_norm – The spectral norm of the inner product. If
None
, it will be computed witheigs
.
- Returns:
Q – The orthonormalized
VectorArray
.R – The upper-triangular/trapezoidal matrix (if
compute_R
isTrue
).