scipy.linalg.pinv#

scipy.linalg.pinv(a, atol=None, rtol=None, return_rank=False, check_finite=True, cond=None, rcond=None)[source]#

Compute the (Moore-Penrose) pseudo-inverse of a matrix.

Calculate a generalized inverse of a matrix using its singular-value decomposition U @ S @ V in the economy mode and picking up only the columns/rows that are associated with significant singular values.

If s is the maximum singular value of a, then the significance cut-off value is determined by atol + rtol * s. Any singular value below this value is assumed insignificant.

Parameters:
a(M, N) array_like

Matrix to be pseudo-inverted.

atolfloat, optional

Absolute threshold term, default value is 0.

New in version 1.7.0.

rtolfloat, optional

Relative threshold term, default value is max(M, N) * eps where eps is the machine precision value of the datatype of a.

New in version 1.7.0.

return_rankbool, optional

If True, return the effective rank of the matrix.

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.

cond, rcondfloat, optional

In older versions, these values were meant to be used as atol with rtol=0. If both were given rcond overwrote cond and hence the code was not correct. Thus using these are strongly discouraged and the tolerances above are recommended instead. In fact, if provided, atol, rtol takes precedence over these keywords.

Changed in version 1.7.0: Deprecated in favor of rtol and atol parameters above and will be removed in future versions of SciPy.

Changed in version 1.3.0: Previously the default cutoff value was just eps*f where f was 1e3 for single precision and 1e6 for double precision.

Returns:
B(N, M) ndarray

The pseudo-inverse of matrix a.

rankint

The effective rank of the matrix. Returned if return_rank is True.

Raises:
LinAlgError

If SVD computation does not converge.

See also

pinvh

Moore-Penrose pseudoinverse of a hermititan matrix.

Notes

If A is invertible then the Moore-Penrose pseudoinverse is exactly the inverse of A [1]. If A is not invertible then the Moore-Penrose pseudoinverse computes the x solution to Ax = b such that ||Ax - b|| is minimized [1].

References

[1] (1,2,3)

Penrose, R. (1956). On best approximate solutions of linear matrix equations. Mathematical Proceedings of the Cambridge Philosophical Society, 52(1), 17-19. doi:10.1017/S0305004100030929

Examples

Given an m x n matrix A and an n x m matrix B the four Moore-Penrose conditions are:

  1. ABA = A (B is a generalized inverse of A),

  2. BAB = B (A is a generalized inverse of B),

  3. (AB)* = AB (AB is hermitian),

  4. (BA)* = BA (BA is hermitian) [1].

Here, A* denotes the conjugate transpose. The Moore-Penrose pseudoinverse is a unique B that satisfies all four of these conditions and exists for any A. Note that, unlike the standard matrix inverse, A does not have to be square or have independant columns/rows.

As an example, we can calculate the Moore-Penrose pseudoinverse of a random non-square matrix and verify it satisfies the four conditions.

>>> import numpy as np
>>> from scipy import linalg
>>> rng = np.random.default_rng()
>>> A = rng.standard_normal((9, 6))
>>> B = linalg.pinv(A)
>>> np.allclose(A @ B @ A, A)  # Condition 1
True
>>> np.allclose(B @ A @ B, B)  # Condition 2
True
>>> np.allclose((A @ B).conj().T, A @ B)  # Condition 3
True
>>> np.allclose((B @ A).conj().T, B @ A)  # Condition 4
True