scipy.linalg.

rq#

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

Compute RQ decomposition of a matrix.

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

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:
a(M, N) array_like

Matrix to be decomposed

overwrite_abool, optional

Whether data in a is overwritten (may improve performance)

lworkint, optional

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

mode{‘full’, ‘r’, ‘economic’}, 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).

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:
Rfloat or complex ndarray

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

Qfloat or complex ndarray

Of shape (N, N) or (K, N) for mode='economic'. Not returned if mode='r'.

Raises:
LinAlgError

If decomposition fails.

Notes

This is an interface to the LAPACK routines sgerqf, dgerqf, cgerqf, zgerqf, sorgrq, dorgrq, cungrq and zungrq.

If mode=economic, the shapes of Q and R are (K, N) and (M, K) instead of (N,N) 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((6, 9))
>>> r, q = linalg.rq(a)
>>> np.allclose(a, r @ q)
True
>>> r.shape, q.shape
((6, 9), (9, 9))
>>> r2 = linalg.rq(a, mode='r')
>>> np.allclose(r, r2)
True
>>> r3, q3 = linalg.rq(a, mode='economic')
>>> r3.shape, q3.shape
((6, 6), (6, 9))