Skip to content

HTTPS clone URL

Subversion checkout URL

You can clone with HTTPS or Subversion.

Download ZIP

Loading…

Pinv enhancements #289

Closed
wants to merge 7 commits into from

5 participants

@jakevdp
Collaborator

This is a PR based on some of @vene's work in scikit-learn: see scikit-learn/scikit-learn#1015

There are two enhancements:

  • scipy.linalg.pinv2 is sped up by only computing necessary singular values, and not allocating the entire psigma matrix

  • scipy.linalg.pinvh added: this uses eigh to significantly speed up the pseudo-inverse computation in the symmetric/hermitian case

Here are some simple benchmarks:

old pinv2 vs new pinv2

old version:

In [1]: import numpy as np
In [2]: from scipy.linalg import pinv2
In [3]: np.random.seed(0)
In [4]: X = np.random.random((500, 100))
In [5]: %timeit pinv2(X)
10 loops, best of 3: 172 ms per loop

new version:

...
In [5]: %timeit pinv2(X)
10 loops, best of 3: 38.1 ms per loop

pinvh vs new pinv2:

In [1]: import numpy as np
In [2]: from scipy.linalg import pinv2, pinvh
In [3]: np.random.seed(0)
In [4]: X = np.random.random((500, 400))
In [5]: X = np.dot(X, X.T) # make symmetric positive semi-definite
In [6]: %timeit pinv2(X)
1 loops, best of 3: 736 ms per loop
In [7]: %timeit pinvh(X)
1 loops, best of 3: 320 ms per loop
@charris charris commented on the diff
scipy/linalg/basic.py
((8 lines not shown))
if rcond is not None:
cond = rcond
if cond in [None,-1]:
- eps = np.finfo(np.float).eps
- feps = np.finfo(np.single).eps
- _array_precision = {'f': 0, 'd': 1, 'F': 0, 'D': 1}
- cond = {0: feps*1e3, 1: eps*1e6}[_array_precision[t]]
- m, n = a.shape
- cutoff = cond*np.maximum.reduce(s)
- psigma = np.zeros((m, n), t)
- for i in range(len(s)):
- if s[i] > cutoff:
- psigma[i,i] = 1.0/np.conjugate(s[i])
- #XXX: use lapack/blas routines for dot
- return np.transpose(np.conjugate(np.dot(np.dot(u,psigma),vh)))
+ t = u.dtype.char.lower()
+ factor = {'f': 1E3, 'd': 1E6}
@charris Collaborator
charris added a note

I wonder if the updated version in np.matrix_rank would be a better choice for the default condition? See the discussion at http://thread.gmane.org/gmane.comp.python.numeric.general/50396/focus=50912

@jakevdp Collaborator
jakevdp added a note

Interesting discussion: here I just duplicated the previous behavior. Any change would be a potential compatibility issue: before making that decision we should probably discuss it on the mailing list.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
@jakevdp
Collaborator

I just had another thought: we should add a lower keyword to pinvh which is passed to eigh.

scipy/linalg/basic.py
((53 lines not shown))
+
+ Returns
+ -------
+ B : array, shape (N, N)
+
+ Raises
+ ------
+ LinAlgError
+ If eigenvalue does not converge
+
+ Examples
+ --------
+ >>> from numpy import *
+ >>> a = random.randn(9, 6)
+ >>> a = np.dot(a, a.T)
+ >>> B = pinv_symmetric(a)

docstring is not up to date right?

@jakevdp Collaborator
jakevdp added a note

right - good catch

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
@jakevdp
Collaborator

@pv - I think this is ready for merge. @charris's comment about matrix rank would be interesting to explore, but I think it's an enhancement that would affect more than just this PR.

scipy/linalg/basic.py
((70 lines not shown))
+ True
+ >>> allclose(B, dot(B, dot(a, B)))
+ True
+
+ """
+ a = np.asarray_chkfinite(a)
+ s, u = decomp.eigh(a, lower=lower)
+
+ if rcond is not None:
+ cond = rcond
+ if cond in [None, -1]:
+ t = u.dtype.char.lower()
+ factor = {'f': 1E3, 'd': 1E6}
+ cond = factor[t] * np.finfo(t).eps
+
+ # unlike svd case, eigh can lead to negative eigenvalues
@pv Owner
pv added a note

More illuminating comment: "For hermitian matrices, singular values equal abs(eigenvalues)"

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
@pv
Owner
pv commented

@jakevdp: looks good to me. Changing the default cutoff is maybe out of scope for this PR. It seems we don't have tests for single-precision data types?

@charris
Collaborator

It might be worth bringing over some of Matthew's documentation for np.matrix_rank, I think it is quite good. But I agree that changing the current default choices for the cutoff is independent of this PR.

@josef-pkt
Collaborator

a suggestion: can we add a keyword argument to optionally return the rank that was used, or eigen or singular values?

looking briefly through the scikit learn PR
"Yes... though I don't think it's possible to do this quickly using the current pinv or pinv2. We should add a flag that optionally returns the number of singular values used in the computation."

In statsmodels, we have currently the problem that we are doing two svd, one for the pinv and one to get the rank. The only option I saw so far was to copy pinv directly into statsmodels with the additional return.

@jakevdp
Collaborator

@josef-pkt: good suggestion on the flag to return the rank. It's a few lines of code that add a lot of flexibility to the use of the routines.

I think I've addressed all the comments (short of the re-working of singular value cutoff selection). This should be ready to merge.

@pv pv referenced this pull request from a commit
@pv pv ENH: linalg: Merge PR #289 Pinv enhancements
This is a PR based on some of @vene's work in scikit-learn: see
scikit-learn/scikit-learn#1015

There are three enhancements:

- scipy.linalg.pinv2 is sped up by only computing necessary singular
  values, and not allocating the entire psigma matrix

- scipy.linalg.pinvh added: this uses eigh to significantly speed up the
  pseudo-inverse computation in the symmetric/hermitian case

- add return_rank keyword to obtain the effective rank from all pinv*
  routines.
5adeb12
@pv
Owner
pv commented

Thanks, merged.

@pv pv closed this
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Commits on Aug 16, 2012
  1. @jakevdp

    speed up linalg.pinv2

    jakevdp authored
  2. @jakevdp
  3. @jakevdp
  4. @jakevdp

    more helpful name for test

    jakevdp authored
  5. @jakevdp

    add lower keyword to pinvh

    jakevdp authored
Commits on Aug 18, 2012
  1. @jakevdp

    cleanup & typos in tests

    jakevdp authored
Commits on Sep 14, 2012
  1. @jakevdp
This page is out of date. Refresh to see the latest.
View
1  scipy/linalg/__init__.py
@@ -31,6 +31,7 @@
lstsq - Solve a linear least-squares problem
pinv - Pseudo-inverse (Moore-Penrose) using lstsq
pinv2 - Pseudo-inverse using svd
+ pinvh - Pseudo-inverse of hermitian matrix
kron - Kronecker product of two arrays
tril - Construct a lower-triangular matrix from a given matrix
triu - Construct an upper-triangular matrix from a given matrix
View
123 scipy/linalg/basic.py
@@ -2,9 +2,10 @@
# Author: Pearu Peterson, March 2002
#
# w/ additions by Travis Oliphant, March 2002
+# and Jake Vanderplas, August 2012
__all__ = ['solve', 'solve_triangular', 'solveh_banded', 'solve_banded',
- 'inv', 'det', 'lstsq', 'pinv', 'pinv2']
+ 'inv', 'det', 'lstsq', 'pinv', 'pinv2', 'pinvh']
import numpy as np
@@ -13,7 +14,7 @@
from misc import LinAlgError, _datacopied
from scipy.linalg import calc_lwork
from decomp_schur import schur
-import decomp_svd
+import decomp, decomp_svd
# Linear equations
@@ -500,7 +501,7 @@ def lstsq(a, b, cond=None, overwrite_a=False, overwrite_b=False):
return x, resids, rank, s
-def pinv(a, cond=None, rcond=None):
+def pinv(a, cond=None, rcond=None, return_rank=False):
"""
Compute the (Moore-Penrose) pseudo-inverse of a matrix.
@@ -515,11 +516,15 @@ def pinv(a, cond=None, rcond=None):
Cutoff for 'small' singular values in the least-squares solver.
Singular values smaller than ``rcond * largest_singular_value``
are considered zero.
+ return_rank : bool, optional
+ if True, return the effective rank of the matrix
Returns
-------
B : array, shape (N, M)
The pseudo-inverse of matrix `a`.
+ rank : int
+ The effective rank of the matrix. Returned if return_rank == True
Raises
------
@@ -540,10 +545,16 @@ def pinv(a, cond=None, rcond=None):
b = np.identity(a.shape[0], dtype=a.dtype)
if rcond is not None:
cond = rcond
- return lstsq(a, b, cond=cond)[0]
+
+ x, resids, rank, s = lstsq(a, b, cond=cond)
+
+ if return_rank:
+ return x, rank
+ else:
+ return x
-def pinv2(a, cond=None, rcond=None):
+def pinv2(a, cond=None, rcond=None, return_rank=False):
"""
Compute the (Moore-Penrose) pseudo-inverse of a matrix.
@@ -560,11 +571,15 @@ def pinv2(a, cond=None, rcond=None):
Singular values smaller than ``rcond*largest_singular_value``
are considered zero.
If None or -1, suitable machine precision is used.
+ return_rank : bool, optional
+ if True, return the effective rank of the matrix
Returns
-------
B : array, shape (N, M)
The pseudo-inverse of matrix `a`.
+ rank : int
+ The effective rank of the matrix. Returned if return_rank == True
Raises
------
@@ -582,21 +597,91 @@ def pinv2(a, cond=None, rcond=None):
"""
a = np.asarray_chkfinite(a)
- u, s, vh = decomp_svd.svd(a)
- t = u.dtype.char
+ u, s, vh = decomp_svd.svd(a, full_matrices=False)
+
if rcond is not None:
cond = rcond
if cond in [None,-1]:
- eps = np.finfo(np.float).eps
- feps = np.finfo(np.single).eps
- _array_precision = {'f': 0, 'd': 1, 'F': 0, 'D': 1}
- cond = {0: feps*1e3, 1: eps*1e6}[_array_precision[t]]
- m, n = a.shape
- cutoff = cond*np.maximum.reduce(s)
- psigma = np.zeros((m, n), t)
- for i in range(len(s)):
- if s[i] > cutoff:
- psigma[i,i] = 1.0/np.conjugate(s[i])
- #XXX: use lapack/blas routines for dot
- return np.transpose(np.conjugate(np.dot(np.dot(u,psigma),vh)))
+ t = u.dtype.char.lower()
+ factor = {'f': 1E3, 'd': 1E6}
@charris Collaborator
charris added a note

I wonder if the updated version in np.matrix_rank would be a better choice for the default condition? See the discussion at http://thread.gmane.org/gmane.comp.python.numeric.general/50396/focus=50912

@jakevdp Collaborator
jakevdp added a note

Interesting discussion: here I just duplicated the previous behavior. Any change would be a potential compatibility issue: before making that decision we should probably discuss it on the mailing list.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
+ cond = factor[t] * np.finfo(t).eps
+
+ above_cutoff = (s > cond * np.max(s))
+ psigma_diag = np.zeros_like(s)
+ psigma_diag[above_cutoff] = 1.0 / s[above_cutoff]
+
+ B = np.transpose(np.conjugate(np.dot(u * psigma_diag, vh)))
+
+ if return_rank:
+ return B, np.sum(above_cutoff)
+ else:
+ return B
+
+
+def pinvh(a, cond=None, rcond=None, lower=True, return_rank=False):
+ """Compute the (Moore-Penrose) pseudo-inverse of a hermetian matrix.
+
+ Calculate a generalized inverse of a symmetric matrix using its
+ eigenvalue decomposition and including all 'large' eigenvalues.
+
+ Parameters
+ ----------
+ a : array, shape (N, N)
+ Real symmetric or complex hermetian matrix to be pseudo-inverted
+ cond, rcond : float or None
+ Cutoff for 'small' eigenvalues.
+ Singular values smaller than rcond * largest_eigenvalue are considered
+ zero.
+
+ If None or -1, suitable machine precision is used.
+ lower : boolean
+ Whether the pertinent array data is taken from the lower or upper
+ triangle of a. (Default: lower)
+ return_rank : bool, optional
+ if True, return the effective rank of the matrix
+
+ Returns
+ -------
+ B : array, shape (N, N)
+ The pseudo-inverse of matrix `a`.
+ rank : int
+ The effective rank of the matrix. Returned if return_rank == True
+
+ Raises
+ ------
+ LinAlgError
+ If eigenvalue does not converge
+
+ Examples
+ --------
+ >>> from numpy import *
+ >>> a = random.randn(9, 6)
+ >>> a = np.dot(a, a.T)
+ >>> B = pinvh(a)
+ >>> allclose(a, dot(a, dot(B, a)))
+ True
+ >>> allclose(B, dot(B, dot(a, B)))
+ True
+ """
+ a = np.asarray_chkfinite(a)
+ s, u = decomp.eigh(a, lower=lower)
+
+ if rcond is not None:
+ cond = rcond
+ if cond in [None, -1]:
+ t = u.dtype.char.lower()
+ factor = {'f': 1E3, 'd': 1E6}
+ cond = factor[t] * np.finfo(t).eps
+
+ # For hermitian matrices, singular values equal abs(eigenvalues)
+ above_cutoff = (abs(s) > cond * np.max(abs(s)))
+ psigma_diag = np.zeros_like(s)
+ psigma_diag[above_cutoff] = 1.0 / s[above_cutoff]
+
+ B = np.dot(u * psigma_diag, np.conjugate(u).T)
+
+ if return_rank:
+ return B, np.sum(above_cutoff)
+ else:
+ return B
View
56 scipy/linalg/tests/test_basic.py
@@ -28,7 +28,7 @@
assert_equal, assert_almost_equal, assert_array_almost_equal, assert_, \
assert_allclose
-from scipy.linalg import solve, inv, det, lstsq, pinv, pinv2, norm,\
+from scipy.linalg import solve, inv, det, lstsq, pinv, pinv2, pinvh, norm,\
solve_banded, solveh_banded, solve_triangular
from scipy.linalg._testutils import assert_no_overwrite
@@ -534,30 +534,66 @@ def test_random_complex_overdet(self):
class TestPinv(TestCase):
- def test_simple(self):
- a=array([[1,2,3],[4,5,6.],[7,8,10]])
+ def test_simple_real(self):
+ a = array([[1, 2, 3], [4, 5, 6], [7, 8, 10]], dtype=float)
a_pinv = pinv(a)
- assert_array_almost_equal(dot(a,a_pinv),[[1,0,0],[0,1,0],[0,0,1]])
+ assert_array_almost_equal(dot(a,a_pinv), np.eye(3))
a_pinv = pinv2(a)
- assert_array_almost_equal(dot(a,a_pinv),[[1,0,0],[0,1,0],[0,0,1]])
- def test_simple_0det(self):
- a=array([[1,2,3],[4,5,6.],[7,8,9]])
+ assert_array_almost_equal(dot(a,a_pinv), np.eye(3))
+
+ def test_simple_complex(self):
+ a = (array([[1, 2, 3], [4, 5, 6], [7, 8, 10]], dtype=float)
+ + 1j * array([[10, 8, 7], [6, 5, 4], [3, 2, 1]], dtype=float))
+ a_pinv = pinv(a)
+ assert_array_almost_equal(dot(a, a_pinv), np.eye(3))
+ a_pinv = pinv2(a)
+ assert_array_almost_equal(dot(a, a_pinv), np.eye(3))
+
+ def test_simple_singular(self):
+ a = array([[1, 2, 3], [4, 5, 6], [7, 8, 9]], dtype=float)
a_pinv = pinv(a)
a_pinv2 = pinv2(a)
assert_array_almost_equal(a_pinv,a_pinv2)
def test_simple_cols(self):
- a=array([[1,2,3],[4,5,6.]])
+ a = array([[1, 2, 3], [4, 5, 6]], dtype=float)
a_pinv = pinv(a)
a_pinv2 = pinv2(a)
assert_array_almost_equal(a_pinv,a_pinv2)
def test_simple_rows(self):
- a=array([[1,2],[3,4],[5,6]])
+ a = array([[1, 2], [3, 4], [5, 6]], dtype=float)
a_pinv = pinv(a)
a_pinv2 = pinv2(a)
assert_array_almost_equal(a_pinv,a_pinv2)
+
+class TestPinvSymmetric(TestCase):
+
+ def test_simple_real(self):
+ a = array([[1, 2, 3], [4, 5, 6], [7, 8, 10]], dtype=float)
+ a = np.dot(a, a.T)
+ a_pinv = pinvh(a)
+ assert_array_almost_equal(np.dot(a, a_pinv), np.eye(3))
+
+ def test_nonpositive(self):
+ a = array([[1, 2, 3], [4, 5, 6], [7, 8, 9]], dtype=float)
+ a = np.dot(a, a.T)
+ u, s, vt = np.linalg.svd(a)
+ s[0] *= -1
+ a = np.dot(u * s, vt) # a is now symmetric non-positive and singular
+ a_pinv = pinv2(a)
+ a_pinvh = pinvh(a)
+ assert_array_almost_equal(a_pinv, a_pinvh)
+
+ def test_simple_complex(self):
+ a = (array([[1, 2, 3], [4, 5, 6], [7, 8, 10]], dtype=float)
+ + 1j * array([[10, 8, 7], [6, 5, 4], [3, 2, 1]], dtype=float))
+ a = np.dot(a, a.conj().T)
+ a_pinv = pinvh(a)
+ assert_array_almost_equal(np.dot(a, a_pinv), np.eye(3))
+
+
class TestNorm(object):
def test_types(self):
@@ -608,6 +644,8 @@ def test_pinv(self):
assert_no_overwrite(pinv, [(3,3)])
def test_pinv2(self):
assert_no_overwrite(pinv2, [(3,3)])
+ def test_pinvh(self):
+ assert_no_overwrite(pinvh, [(3,3)])
if __name__ == "__main__":
run_module_suite()
Something went wrong with that request. Please try again.