scipy.linalg.pinv#
- scipy.linalg.pinv(a, *, atol=None, rtol=None, return_rank=False, check_finite=True, cond=<object object>, rcond=<object object>)[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 ofa
, then the significance cut-off value is determined byatol + 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
whereeps
is the machine precision value of the datatype ofa
.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
withrtol=0
. If both were givenrcond
overwrotecond
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.Deprecated since version 1.7.0: Deprecated in favor of
rtol
andatol
parameters above and will be removed in SciPy 1.14.0.Changed in version 1.3.0: Previously the default cutoff value was just
eps*f
wheref
was1e3
for single precision and1e6
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 hermitian matrix.
Notes
If
A
is invertible then the Moore-Penrose pseudoinverse is exactly the inverse ofA
[1]. IfA
is not invertible then the Moore-Penrose pseudoinverse computes thex
solution toAx = b
such that||Ax - b||
is minimized [1].References
Examples
Given an
m x n
matrixA
and ann x m
matrixB
the four Moore-Penrose conditions are:ABA = A
(B
is a generalized inverse ofA
),BAB = B
(A
is a generalized inverse ofB
),(AB)* = AB
(AB
is hermitian),(BA)* = BA
(BA
is hermitian) [1].
Here,
A*
denotes the conjugate transpose. The Moore-Penrose pseudoinverse is a uniqueB
that satisfies all four of these conditions and exists for anyA
. Note that, unlike the standard matrix inverse,A
does not have to be a square matrix or have linearly independent 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