qr_insert#
- scipy.linalg.qr_insert(Q, R, u, k, which='row', rcond=None, overwrite_qru=False, check_finite=True)#
QR update on row or column insertions
If
A = Q R
is the QR factorization ofA
, return the QR factorization ofA
where rows or columns have been inserted starting at row or columnk
.The documentation is written assuming array arguments are of specified “core” shapes. However, array argument(s) of this function may have additional “batch” dimensions prepended to the core shape. In this case, the array is treated as a batch of lower-dimensional slices; see Batched Linear Operations for details.
- Parameters:
- Q(M, M) array_like
Unitary/orthogonal matrix from the QR decomposition of A.
- R(M, N) array_like
Upper triangular matrix from the QR decomposition of A.
- u(N,), (p, N), (M,), or (M, p) array_like
Rows or columns to insert
- kint
Index before which u is to be inserted.
- which: {‘row’, ‘col’}, optional
Determines if rows or columns will be inserted, defaults to ‘row’
- rcondfloat
Lower bound on the reciprocal condition number of
Q
augmented withu/||u||
Only used when updating economic mode (thin, (M,N) (N,N)) decompositions. If None, machine precision is used. Defaults to None.- overwrite_qrubool, optional
If True, consume Q, R, and u, if possible, while performing the update, otherwise make copies as necessary. Defaults to False.
- check_finitebool, optional
Whether to check that the input matrices contain only finite numbers. Disabling may give a performance gain, but may result in problems (crashes, non-termination) if the inputs do contain infinities or NaNs. Default is True.
- Returns:
- Q1ndarray
Updated unitary/orthogonal factor
- R1ndarray
Updated upper triangular factor
- Raises:
- LinAlgError
If updating a (M,N) (N,N) factorization and the reciprocal condition number of Q augmented with
u/||u||
is smaller than rcond.
See also
Notes
This routine does not guarantee that the diagonal entries of
R1
are positive.Added in version 0.16.0.
References
[1]Golub, G. H. & Van Loan, C. F. Matrix Computations, 3rd Ed. (Johns Hopkins University Press, 1996).
[2]Daniel, J. W., Gragg, W. B., Kaufman, L. & Stewart, G. W. Reorthogonalization and stable algorithms for updating the Gram-Schmidt QR factorization. Math. Comput. 30, 772-795 (1976).
[3]Reichel, L. & Gragg, W. B. Algorithm 686: FORTRAN Subroutines for Updating the QR Decomposition. ACM Trans. Math. Softw. 16, 369-377 (1990).
Examples
>>> import numpy as np >>> from scipy import linalg >>> a = np.array([[ 3., -2., -2.], ... [ 6., -7., 4.], ... [ 7., 8., -6.]]) >>> q, r = linalg.qr(a)
Given this QR decomposition, update q and r when 2 rows are inserted.
>>> u = np.array([[ 6., -9., -3.], ... [ -3., 10., 1.]]) >>> q1, r1 = linalg.qr_insert(q, r, u, 2, 'row') >>> q1 array([[-0.25445668, 0.02246245, 0.18146236, -0.72798806, 0.60979671], # may vary (signs) [-0.50891336, 0.23226178, -0.82836478, -0.02837033, -0.00828114], [-0.50891336, 0.35715302, 0.38937158, 0.58110733, 0.35235345], [ 0.25445668, -0.52202743, -0.32165498, 0.36263239, 0.65404509], [-0.59373225, -0.73856549, 0.16065817, -0.0063658 , -0.27595554]]) >>> r1 array([[-11.78982612, 6.44623587, 3.81685018], # may vary (signs) [ 0. , -16.01393278, 3.72202865], [ 0. , 0. , -6.13010256], [ 0. , 0. , 0. ], [ 0. , 0. , 0. ]])
The update is equivalent, but faster than the following.
>>> a1 = np.insert(a, 2, u, 0) >>> a1 array([[ 3., -2., -2.], [ 6., -7., 4.], [ 6., -9., -3.], [ -3., 10., 1.], [ 7., 8., -6.]]) >>> q_direct, r_direct = linalg.qr(a1)
Check that we have equivalent results:
>>> np.dot(q1, r1) array([[ 3., -2., -2.], [ 6., -7., 4.], [ 6., -9., -3.], [ -3., 10., 1.], [ 7., 8., -6.]])
>>> np.allclose(np.dot(q1, r1), a1) True
And the updated Q is still unitary:
>>> np.allclose(np.dot(q1.T, q1), np.eye(5)) True