Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

ENH: Vectorising np.linalg.qr #19151

Merged
merged 26 commits into from Jul 14, 2021
Merged
Show file tree
Hide file tree
Changes from 25 commits
Commits
Show all changes
26 commits
Select commit Hold shift + click to select a range
c25c09b
ENH: r, raw, economic now accept stacked matrices
czgdp1807 Jun 3, 2021
d5cbb5e
WIP
czgdp1807 Jun 4, 2021
a3b58ba
ready for formal addition of tests and documentation update
czgdp1807 Jun 4, 2021
885ba91
moving intialization outside the loop
czgdp1807 Jun 5, 2021
c33cf14
existing tests passed
czgdp1807 Jun 5, 2021
b6cd5b2
tests for stacked inputs added
czgdp1807 Jun 5, 2021
4da4e5e
documentation updated
czgdp1807 Jun 5, 2021
89ff2af
linting issue fixed
czgdp1807 Jun 5, 2021
4628497
removed trailing white spaces
czgdp1807 Jun 5, 2021
9418a17
resolved linting issues
czgdp1807 Jun 5, 2021
0d490d6
resolved linting issues
czgdp1807 Jun 5, 2021
3b64375
documentation fix
czgdp1807 Jun 5, 2021
cd66940
removed debug prints and addressed review
czgdp1807 Jun 7, 2021
3e7baee
Added release notes entry
czgdp1807 Jun 7, 2021
c2dc5f7
Addressed reviews
czgdp1807 Jun 8, 2021
7529874
fixed title line in release notes
czgdp1807 Jun 8, 2021
38fcfc2
Take work_count as it is
czgdp1807 Jun 9, 2021
03bfafc
shifted qr code in umath_linalg.c.src above
czgdp1807 Jun 9, 2021
3f1c88c
updated method for finding optimal work_count
czgdp1807 Jun 9, 2021
e9df9fd
factored out common init parts for sorgqr, dorgqr
czgdp1807 Jun 22, 2021
b5d469e
refactoring complete
czgdp1807 Jun 22, 2021
b498d65
Addressed testing reviews
czgdp1807 Jun 24, 2021
e88849d
addressed C reviews
czgdp1807 Jun 24, 2021
c1f500b
fixed linting issues
czgdp1807 Jun 24, 2021
6ab34ef
removed redundant functions
czgdp1807 Jul 9, 2021
6e405d5
removed redudancy in code
czgdp1807 Jul 10, 2021
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
6 changes: 6 additions & 0 deletions doc/release/upcoming_changes/19151.improvement.rst
@@ -0,0 +1,6 @@
`numpy.linalg.qr` accepts stacked matrices as inputs
----------------------------------------------------

`numpy.linalg.qr` is able to produce results for stacked matrices as inputs.
Moreover, the implementation of QR decomposition has been shifted to C
from Python.
130 changes: 66 additions & 64 deletions numpy/linalg/linalg.py
Expand Up @@ -99,6 +99,10 @@ def _raise_linalgerror_svd_nonconvergence(err, flag):
def _raise_linalgerror_lstsq(err, flag):
raise LinAlgError("SVD did not converge in Linear Least Squares")

def _raise_linalgerror_qr(err, flag):
raise LinAlgError("Incorrect argument found while performing "
"QR factorization")

def get_linalg_error_extobj(callback):
extobj = list(_linalg_error_extobj) # make a copy
extobj[2] = callback
Expand Down Expand Up @@ -780,15 +784,16 @@ def qr(a, mode='reduced'):

Parameters
----------
a : array_like, shape (M, N)
Matrix to be factored.
a : array_like, shape (..., M, N)
An array-like object with the dimensionality of at least 2.
mode : {'reduced', 'complete', 'r', 'raw'}, optional
If K = min(M, N), then

* 'reduced' : returns q, r with dimensions (M, K), (K, N) (default)
* 'complete' : returns q, r with dimensions (M, M), (M, N)
* 'r' : returns r only with dimensions (K, N)
* 'raw' : returns h, tau with dimensions (N, M), (K,)
* 'reduced' : returns q, r with dimensions
(..., M, K), (..., K, N) (default)
* 'complete' : returns q, r with dimensions (..., M, M), (..., M, N)
* 'r' : returns r only with dimensions (..., K, N)
* 'raw' : returns h, tau with dimensions (..., N, M), (..., K,)

The options 'reduced', 'complete, and 'raw' are new in numpy 1.8,
see the notes for more information. The default is 'reduced', and to
Expand All @@ -807,9 +812,13 @@ def qr(a, mode='reduced'):
A matrix with orthonormal columns. When mode = 'complete' the
result is an orthogonal/unitary matrix depending on whether or not
a is real/complex. The determinant may be either +/- 1 in that
case.
case. In case the number of dimensions in the input array is
greater than 2 then a stack of the matrices with above properties
is returned.
r : ndarray of float or complex, optional
The upper-triangular matrix.
The upper-triangular matrix or a stack of upper-triangular
matrices if the number of dimensions in the input array is greater
than 2.
(h, tau) : ndarrays of np.double or np.cdouble, optional
The array h contains the Householder reflectors that generate q
along with r. The tau array contains scaling factors for the
Expand Down Expand Up @@ -857,6 +866,14 @@ def qr(a, mode='reduced'):
>>> r2 = np.linalg.qr(a, mode='r')
>>> np.allclose(r, r2) # mode='r' returns the same r as mode='full'
True
>>> a = np.random.normal(size=(3, 2, 2)) # Stack of 2 x 2 matrices as input
>>> q, r = np.linalg.qr(a)
>>> q.shape
(3, 2, 2)
>>> r.shape
(3, 2, 2)
>>> np.allclose(a, np.matmul(q, r))
True

Example illustrating a common use of `qr`: solving of least squares
problems
Expand Down Expand Up @@ -904,83 +921,68 @@ def qr(a, mode='reduced'):
raise ValueError(f"Unrecognized mode '{mode}'")

a, wrap = _makearray(a)
_assert_2d(a)
m, n = a.shape
_assert_stacked_2d(a)
m, n = a.shape[-2:]
t, result_t = _commonType(a)
a = _fastCopyAndTranspose(t, a)
a = a.astype(t, copy=True)
a = _to_native_byte_order(a)
mn = min(m, n)
tau = zeros((mn,), t)

if isComplexType(t):
lapack_routine = lapack_lite.zgeqrf
routine_name = 'zgeqrf'
if m <= n:
gufunc = _umath_linalg.qr_r_raw_m
else:
lapack_routine = lapack_lite.dgeqrf
routine_name = 'dgeqrf'

# calculate optimal size of work data 'work'
lwork = 1
work = zeros((lwork,), t)
results = lapack_routine(m, n, a, max(1, m), tau, work, -1, 0)
if results['info'] != 0:
raise LinAlgError('%s returns %d' % (routine_name, results['info']))

# do qr decomposition
lwork = max(1, n, int(abs(work[0])))
work = zeros((lwork,), t)
results = lapack_routine(m, n, a, max(1, m), tau, work, lwork, 0)
if results['info'] != 0:
raise LinAlgError('%s returns %d' % (routine_name, results['info']))
gufunc = _umath_linalg.qr_r_raw_n

signature = 'D->D' if isComplexType(t) else 'd->d'
mattip marked this conversation as resolved.
Show resolved Hide resolved
extobj = get_linalg_error_extobj(_raise_linalgerror_qr)
tau = gufunc(a, signature=signature, extobj=extobj)

# handle modes that don't return q
if mode == 'r':
r = _fastCopyAndTranspose(result_t, a[:, :mn])
return wrap(triu(r))
r = triu(a[..., :mn, :])
if t != result_t:
r = r.astype(result_t, copy=False)
pearu marked this conversation as resolved.
Show resolved Hide resolved
return wrap(r)

if mode == 'raw':
return a, tau
q = transpose(a)
if t != result_t:
pearu marked this conversation as resolved.
Show resolved Hide resolved
q = q.astype(result_t, copy=False)
tau = tau.astype(result_t, copy=False)
return wrap(q), tau

if mode == 'economic':
if t != result_t :
pearu marked this conversation as resolved.
Show resolved Hide resolved
a = a.astype(result_t, copy=False)
return wrap(a.T)
return wrap(a)
pearu marked this conversation as resolved.
Show resolved Hide resolved

# generate q from a
# mc is the number of columns in the resulting q
# matrix. If the mode is complete then it is
# same as number of rows, and if the mode is reduced,
# then it is the minimum of number of rows and columns.
if mode == 'complete' and m > n:
mc = m
q = empty((m, m), t)
if m <= n:
gufunc = _umath_linalg.qr_complete_m
else:
gufunc = _umath_linalg.qr_complete_n
pearu marked this conversation as resolved.
Show resolved Hide resolved
else:
mc = mn
q = empty((n, m), t)
q[:n] = a

if isComplexType(t):
lapack_routine = lapack_lite.zungqr
routine_name = 'zungqr'
else:
lapack_routine = lapack_lite.dorgqr
routine_name = 'dorgqr'

# determine optimal lwork
lwork = 1
work = zeros((lwork,), t)
results = lapack_routine(m, mc, mn, q, max(1, m), tau, work, -1, 0)
if results['info'] != 0:
raise LinAlgError('%s returns %d' % (routine_name, results['info']))

# compute q
lwork = max(1, n, int(abs(work[0])))
work = zeros((lwork,), t)
results = lapack_routine(m, mc, mn, q, max(1, m), tau, work, lwork, 0)
if results['info'] != 0:
raise LinAlgError('%s returns %d' % (routine_name, results['info']))
if m <= n:
gufunc = _umath_linalg.qr_reduced_m
else:
gufunc = _umath_linalg.qr_reduced_n

q = _fastCopyAndTranspose(result_t, q[:mc])
r = _fastCopyAndTranspose(result_t, a[:, :mc])
signature = 'DD->D' if isComplexType(t) else 'dd->d'
extobj = get_linalg_error_extobj(_raise_linalgerror_qr)
q = gufunc(a, tau, signature=signature, extobj=extobj)
r = triu(a[..., :mc, :])
czgdp1807 marked this conversation as resolved.
Show resolved Hide resolved

return wrap(q), wrap(triu(r))
if t != result_t:
pearu marked this conversation as resolved.
Show resolved Hide resolved
q = q.astype(result_t, copy=False)
r = r.astype(result_t, copy=False)

return wrap(q), wrap(r)

# Eigenvalues

Expand Down Expand Up @@ -2177,7 +2179,7 @@ def lstsq(a, b, rcond="warn"):
equal to, or greater than its number of linearly independent columns).
If `a` is square and of full rank, then `x` (but for round-off error)
is the "exact" solution of the equation. Else, `x` minimizes the
Euclidean 2-norm :math:`||b - ax||`. If there are multiple minimizing
Euclidean 2-norm :math:`||b - ax||`. If there are multiple minimizing
solutions, the one with the smallest 2-norm :math:`||x||` is returned.

Parameters
Expand Down
62 changes: 62 additions & 0 deletions numpy/linalg/tests/test_linalg.py
@@ -1,6 +1,7 @@
""" Test functions for linalg module

"""
from numpy.core.fromnumeric import shape
import os
import sys
import itertools
Expand All @@ -11,6 +12,7 @@

import numpy as np
from numpy import array, single, double, csingle, cdouble, dot, identity, matmul
from numpy.core import swapaxes
from numpy import multiply, atleast_2d, inf, asarray
from numpy import linalg
from numpy.linalg import matrix_power, norm, matrix_rank, multi_dot, LinAlgError
Expand Down Expand Up @@ -1710,6 +1712,66 @@ def test_mode_all_but_economic(self):
self.check_qr(m2)
self.check_qr(m2.T)

def check_qr_stacked(self, a):
# This test expects the argument `a` to be an ndarray or
# a subclass of an ndarray of inexact type.
a_type = type(a)
a_dtype = a.dtype
m, n = a.shape[-2:]
k = min(m, n)

# mode == 'complete'
q, r = linalg.qr(a, mode='complete')
assert_(q.dtype == a_dtype)
assert_(r.dtype == a_dtype)
assert_(isinstance(q, a_type))
assert_(isinstance(r, a_type))
assert_(q.shape[-2:] == (m, m))
assert_(r.shape[-2:] == (m, n))
assert_almost_equal(matmul(q, r), a)
I_mat = np.identity(q.shape[-1])
stack_I_mat = np.broadcast_to(I_mat,
q.shape[:-2] + (q.shape[-1],)*2)
assert_almost_equal(matmul(swapaxes(q, -1, -2).conj(), q), stack_I_mat)
assert_almost_equal(np.triu(r[..., :, :]), r)

# mode == 'reduced'
q1, r1 = linalg.qr(a, mode='reduced')
assert_(q1.dtype == a_dtype)
assert_(r1.dtype == a_dtype)
assert_(isinstance(q1, a_type))
assert_(isinstance(r1, a_type))
assert_(q1.shape[-2:] == (m, k))
assert_(r1.shape[-2:] == (k, n))
assert_almost_equal(matmul(q1, r1), a)
I_mat = np.identity(q1.shape[-1])
stack_I_mat = np.broadcast_to(I_mat,
q1.shape[:-2] + (q1.shape[-1],)*2)
assert_almost_equal(matmul(swapaxes(q1, -1, -2).conj(), q1),
stack_I_mat)
assert_almost_equal(np.triu(r1[..., :, :]), r1)

# mode == 'r'
r2 = linalg.qr(a, mode='r')
assert_(r2.dtype == a_dtype)
assert_(isinstance(r2, a_type))
assert_almost_equal(r2, r1)

@pytest.mark.parametrize("size", [
(3, 4), (4, 3), (4, 4),
(3, 0), (0, 3)])
@pytest.mark.parametrize("outer_size", [
(2, 2), (2,), (2, 3, 4)])
@pytest.mark.parametrize("dt", [
np.single, np.double,
np.csingle, np.cdouble])
def test_stacked_inputs(self, outer_size, size, dt):

A = np.random.normal(size=outer_size + size).astype(dt)
B = np.random.normal(size=outer_size + size).astype(dt)
self.check_qr_stacked(A)
self.check_qr_stacked(A + 1.j*B)


class TestCholesky:
# TODO: are there no other tests for cholesky?
Expand Down