scipy.linalg.qr#

scipy.linalg.qr(a, overwrite_a=False, lwork=None, mode='full', pivoting=False, check_finite=True)[source]#

Compute QR decomposition of a matrix.

Calculate the decomposition A = Q R where Q is unitary/orthogonal and R upper triangular.

Parameters:
a(M, N) array_like

Matrix to be decomposed

overwrite_abool, optional

Whether data in a is overwritten (may improve performance if overwrite_a is set to True by reusing the existing input data structure rather than creating a new one.)

lworkint, optional

Work array size, lwork >= a.shape[1]. If None or -1, an optimal size is computed.

mode{‘full’, ‘r’, ‘economic’, ‘raw’}, optional

Determines what information is to be returned: either both Q and R (‘full’, default), only R (‘r’) or both Q and R but computed in economy-size (‘economic’, see Notes). The final option ‘raw’ (added in SciPy 0.11) makes the function return two matrices (Q, TAU) in the internal format used by LAPACK.

pivotingbool, optional

Whether or not factorization should include pivoting for rank-revealing qr decomposition. If pivoting, compute the decomposition A[:, P] = Q @ R as above, but where P is chosen such that the diagonal of R is non-increasing.

check_finitebool, optional

Whether to check that the input matrix contains 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.

Returns:
Qfloat or complex ndarray

Of shape (M, M), or (M, K) for mode='economic'. Not returned if mode='r'. Replaced by tuple (Q, TAU) if mode='raw'.

Rfloat or complex ndarray

Of shape (M, N), or (K, N) for mode in ['economic', 'raw']. K = min(M, N).

Pint ndarray

Of shape (N,) for pivoting=True. Not returned if pivoting=False.

Raises:
LinAlgError

Raised if decomposition fails

Notes

This is an interface to the LAPACK routines dgeqrf, zgeqrf, dorgqr, zungqr, dgeqp3, and zgeqp3.

If mode=economic, the shapes of Q and R are (M, K) and (K, N) instead of (M,M) and (M,N), with K=min(M,N).

Examples

>>> import numpy as np
>>> from scipy import linalg
>>> rng = np.random.default_rng()
>>> a = rng.standard_normal((9, 6))
>>> q, r = linalg.qr(a)
>>> np.allclose(a, np.dot(q, r))
True
>>> q.shape, r.shape
((9, 9), (9, 6))
>>> r2 = linalg.qr(a, mode='r')
>>> np.allclose(r, r2)
True
>>> q3, r3 = linalg.qr(a, mode='economic')
>>> q3.shape, r3.shape
((9, 6), (6, 6))
>>> q4, r4, p4 = linalg.qr(a, pivoting=True)
>>> d = np.abs(np.diag(r4))
>>> np.all(d[1:] <= d[:-1])
True
>>> np.allclose(a[:, p4], np.dot(q4, r4))
True
>>> q4.shape, r4.shape, p4.shape
((9, 9), (9, 6), (6,))
>>> q5, r5, p5 = linalg.qr(a, mode='economic', pivoting=True)
>>> q5.shape, r5.shape, p5.shape
((9, 6), (6, 6), (6,))