qr_multiply#
- scipy.linalg.qr_multiply(a, c, mode='right', pivoting=False, conjugate=False, overwrite_a=False, overwrite_c=False)[source]#
Calculate the QR decomposition and multiply Q with a matrix.
Calculate the decomposition
A = Q R
where Q is unitary/orthogonal and R upper triangular. Multiply Q with a vector or a matrix c.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
Input array
- carray_like
Input array to be multiplied by
q
.- mode{‘left’, ‘right’}, optional
Q @ c
is returned if mode is ‘left’,c @ Q
is returned if mode is ‘right’. The shape of c must be appropriate for the matrix multiplications, if mode is ‘left’,min(a.shape) == c.shape[0]
, if mode is ‘right’,a.shape[0] == c.shape[1]
.- pivotingbool, optional
Whether or not factorization should include pivoting for rank-revealing qr decomposition, see the documentation of qr.
- conjugatebool, optional
Whether Q should be complex-conjugated. This might be faster than explicit conjugation.
- overwrite_abool, optional
Whether data in a is overwritten (may improve performance)
- overwrite_cbool, optional
Whether data in c is overwritten (may improve performance). If this is used, c must be big enough to keep the result, i.e.
c.shape[0]
=a.shape[0]
if mode is ‘left’.
- Returns:
- CQndarray
The product of
Q
andc
.- R(K, N), ndarray
R array of the resulting QR factorization where
K = min(M, N)
.- P(N,) ndarray
Integer pivot array. Only returned when
pivoting=True
.
- Raises:
- LinAlgError
Raised if QR decomposition fails.
Notes
This is an interface to the LAPACK routines
?GEQRF
,?ORMQR
,?UNMQR
, and?GEQP3
.Added in version 0.11.0.
Examples
>>> import numpy as np >>> from scipy.linalg import qr_multiply, qr >>> A = np.array([[1, 3, 3], [2, 3, 2], [2, 3, 3], [1, 3, 2]]) >>> qc, r1, piv1 = qr_multiply(A, 2*np.eye(4), pivoting=1) >>> qc array([[-1., 1., -1.], [-1., -1., 1.], [-1., -1., -1.], [-1., 1., 1.]]) >>> r1 array([[-6., -3., -5. ], [ 0., -1., -1.11022302e-16], [ 0., 0., -1. ]]) >>> piv1 array([1, 0, 2], dtype=int32) >>> q2, r2, piv2 = qr(A, mode='economic', pivoting=1) >>> np.allclose(2*q2 - qc, np.zeros((4, 3))) True