scipy.linalg.

null_space#

scipy.linalg.null_space(A, rcond=None, *, overwrite_a=False, check_finite=True, lapack_driver='gesdd')[source]#

Construct an orthonormal basis for the null space of A using SVD

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

rcondfloat, optional

Relative condition number. Singular values s smaller than rcond * max(s) are considered zero. Default: floating point eps * max(M,N).

overwrite_abool, optional

Whether to overwrite a; may improve performance. Default is False.

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.

lapack_driver{‘gesdd’, ‘gesvd’}, optional

Whether to use the more efficient divide-and-conquer approach ('gesdd') or general rectangular approach ('gesvd') to compute the SVD. MATLAB and Octave use the 'gesvd' approach. Default is 'gesdd'.

Returns:
Z(N, K) ndarray

Orthonormal basis for the null space of A. K = dimension of effective null space, as determined by rcond

See also

svd

Singular value decomposition of a matrix

orth

Matrix range

Examples

1-D null space:

>>> import numpy as np
>>> from scipy.linalg import null_space
>>> A = np.array([[1, 1], [1, 1]])
>>> ns = null_space(A)
>>> ns * np.copysign(1, ns[0,0])  # Remove the sign ambiguity of the vector
array([[ 0.70710678],
       [-0.70710678]])

2-D null space:

>>> from numpy.random import default_rng
>>> rng = default_rng()
>>> B = rng.random((3, 5))
>>> Z = null_space(B)
>>> Z.shape
(5, 2)
>>> np.allclose(B.dot(Z), 0)
True

The basis vectors are orthonormal (up to rounding error):

>>> Z.T.dot(Z)
array([[  1.00000000e+00,   6.92087741e-17],
       [  6.92087741e-17,   1.00000000e+00]])