Return an orthonormal basis for the kernel of mat 1.
The routine employs the singular value decomposition so that the columns of the
returned matrix span the null space and are orthonormal with respect to the
standard inner product.
Parameters:
-
mat
(ndarray)
–
Matrix whose null space is sought.
-
tol
(float, default:
1e-08
)
–
Numerical tolerance that distinguishes zero singular values.
Returns:
-
ndarray
–
A matrix whose columns form an orthonormal basis for the null space.
Examples:
Consider the matrix
\[
A = \begin{pmatrix} 1 & 1 & 0 \\ 0 & 0 & 0 \end{pmatrix}.
\]
Its null space is spanned by the vectors \((1,-1,0)\) and \((0,0,1)\).
import numpy as np
from toqito.matrix_ops import null_space
A = np.array([[1, 1, 0], [0, 0, 0]], dtype=float)
null_basis = null_space(A)
print(null_basis)
[[-0.70710678+0.j 0. +0.j]
[ 0.70710678-0.j 0. +0.j]
[ 0. -0.j -1. +0.j]]
References
1 Wikipedia. Kernel (linear algebra). link.
Source code in toqito/matrix_ops/null_space.py
| def null_space(mat: np.ndarray, tol: float = 1e-08) -> np.ndarray:
r"""Return an orthonormal basis for the kernel of ``mat`` [@wikipediakernel].
The routine employs the singular value decomposition so that the columns of the
returned matrix span the null space and are orthonormal with respect to the
standard inner product.
Args:
mat: Matrix whose null space is sought.
tol: Numerical tolerance that distinguishes zero singular values.
Returns:
A matrix whose columns form an orthonormal basis for the null space.
Examples:
Consider the matrix
\[
A = \begin{pmatrix} 1 & 1 & 0 \\ 0 & 0 & 0 \end{pmatrix}.
\]
Its null space is spanned by the vectors \((1,-1,0)\) and \((0,0,1)\).
```python exec="1" source="above" result="text"
import numpy as np
from toqito.matrix_ops import null_space
A = np.array([[1, 1, 0], [0, 0, 0]], dtype=float)
null_basis = null_space(A)
print(null_basis)
```
"""
mat = np.asarray(mat, dtype=np.complex128)
if mat.ndim != 2:
raise ValueError("Input must be a two-dimensional array.")
if mat.size == 0:
return np.zeros((mat.shape[0], 0), dtype=np.complex128)
_, singular_values, vh = np.linalg.svd(mat, full_matrices=True)
rank = np.sum(singular_values > tol)
kernel = vh[rank:].conj().T
if kernel.size == 0:
return np.zeros((mat.shape[1], 0), dtype=np.complex128)
q, _ = np.linalg.qr(kernel)
return q
|