nitorch_fastmath.sugar
Overview
This module contains "syntactic sugar" linear algebra routines. I.e., their behaviour is relatively easily to implement in pure PyTorch, but sometimes not very readable as it involves squeezing, unsqueezing, or shuffling dimensions.
For example, we implement matvec(mat, vec), which in PyTorch must
be written matmul(mat, vec.unsqueeze(-1)).squeeze(-1).
kron2
Kronecker product of two matrices \(\mathbf{A} \otimes \mathbf{B}\)
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
Tensor
|
Left matrix, with shape |
required |
b |
`(..., p, q) tensor`
|
Right matrix, with shape |
required |
Returns:
| Name | Type | Description |
|---|---|---|
ab |
`(..., p*m, q*n) tensor`
|
Kronecker product, with shape Note: |
References
- Kronecker product, Wikipedia.
lmdiv
Left matrix division \(\mathbf{A}^{-1} \times \mathbf{B}\).
Note
if m != n, the Moore-Penrose pseudoinverse is always used.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
`(..., m, n) tensor`
|
Left input ("the system"), with shape |
required |
b |
`(..., m, k) tensor`
|
Right input ("the point"), with shape |
required |
method |
`{'lu', 'chol', 'svd', 'pinv'}`
|
Inversion method:
|
'lu'
|
rcond |
`float`
|
Cutoff for small singular values when |
1e-15
|
out |
`tensor`
|
Output tensor. |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
x |
`(..., n, k) tensor`
|
Solution of the linear system, with shape |
References
- LU decomposition, Wikipedia.
- Cholesky decomposition, Wikipedia.
- Singular value decomposition, Wikipedia.
- Moore-Penrose inverse, Wikipedia.
rmdiv
Right matrix division \(\mathbf{A} \times \mathbf{B}^{-1}\).
Note
if m != n, the Moore-Penrose pseudoinverse is always used.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
`(..., k, m) tensor`
|
Left input ("the point"), with shape |
required |
b |
`(..., n, m) tensor`
|
Right input ("the system"), with shape |
required |
method |
`{'lu', 'chol', 'svd', 'pinv'}`
|
Inversion method:
|
'lu'
|
rcond |
`float`
|
Cutoff for small singular values when |
1e-15
|
out |
`tensor`
|
Output tensor |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
x |
`(..., k, m) tensor`
|
Solution of the linear system, with shape |
References
- LU decomposition, Wikipedia.
- Cholesky decomposition, Wikipedia.
- Singular value decomposition, Wikipedia.
- Moore-Penrose inverse, Wikipedia.
inv
Matrix inversion \(\mathbf{A}^{-1}\).
Note
if m != n, the Moore-Penrose pseudoinverse is always used.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
`(..., m, n) tensor`
|
Input matrix, with shape |
required |
method |
`{'lu', 'chol', 'svd', 'pinv'}`
|
Inversion method:
|
'lu'
|
rcond |
`float`
|
Cutoff for small singular values when |
1e-15
|
out |
`tensor`
|
Output tensor). |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
x |
`(..., n, m) tensor`
|
Inverse matrix, with shape |
References
- LU decomposition, Wikipedia.
- Cholesky decomposition, Wikipedia.
- Singular value decomposition, Wikipedia.
- Moore-Penrose inverse, Wikipedia.
matvec
Matrix-vector product \(\mathbf{A} \times \mathbf{b}\) (supports broadcasting) .
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
mat |
`(..., m, n) tensor`
|
Input matrix, with shape |
required |
vec |
`(..., n) tensor`
|
Input vector, with shape |
required |
out |
`(..., n) tensor`
|
Placeholder for the output tensor, with shape |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
mv |
`(..., m) tensor`
|
Matrix vector product of the inputs, with shape |
solvevec
Left matrix-vector division \(\mathbf{A}^{-1} \times \mathbf{b}\).
Note
if m != n, the Moore-Penrose pseudoinverse is always used.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
mat |
`(..., m, n) tensor`
|
Left input ("the system"), with shape |
required |
vec |
`(..., m) tensor`
|
Right input ("the point"), with shape |
required |
method |
`{'lu', 'chol', 'svd', 'pinv'}`
|
Inversion method:
|
'lu'
|
rcond |
`float`
|
Cutoff for small singular values when |
1e-15
|
out |
`tensor`
|
Output tensor. |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
x |
`(..., n) tensor`
|
Solution of the linear system, with shape |
References
- LU decomposition, Wikipedia.
- Cholesky decomposition, Wikipedia.
- Singular value decomposition, Wikipedia.
- Moore-Penrose inverse, Wikipedia.
outer
Outer product of two (batched) tensors \(\mathbf{a} \otimes \mathbf{b} = \mathbf{a} \times \mathbf{b}^\mathrm{H}\).
Warning
When \(\mathbf{a}\) and \(\mathbf{b}\) are complex, this function returns the complex outer product \(\mathbf{a} \otimes \mathbf{b} = \mathbf{a}\mathbf{b}^\mathrm{H}\).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
`(..., n) tensor`
|
Left vector, with shape |
required |
b |
`(..., m) tensor`
|
Right vector, with shape |
required |
out |
`(..., n, m) tensor`
|
Output placeholder, with shape |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
out |
`(..., n, m) tensor`
|
Outer product, with shape |
References
- Outer product, Wikipedia.
trace
Compute the (batched) trace of a matrix \(\operatorname{tr}\left(\mathbf{A}\right)\).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
`(..., m, m) tensor`
|
|
required |
Returns:
| Name | Type | Description |
|---|---|---|
t |
`(..., [1, 1]) tensor`
|
|
References
- Trace, Wikipedia.
dot
(Batched) dot product of two vectors \(\mathbf{a}^\mathrm{H}\mathbf{b}\).
Warning
When \(\mathbf{a}\) and \(\mathbf{b}\) are complex, this function returns the complex dot product \((\mathbf{a}, \mathbf{b}) = \mathbf{a}^\mathrm{H}\mathbf{b}\). This dot product is linear in the second term and antilinear in the first term.
This behaviour differs from torch.dot and torch.inner, which
both compute \(\mathbf{a}^\mathrm{T}\mathbf{b}\) when the vectors are complex.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
`(..., n) tensor`
|
Left vector, with shape |
required |
b |
`(..., m) tensor`
|
Right vector, with shape |
required |
keepdim |
`bool`
|
Keep reduced dimension |
False
|
out |
`tensor`
|
Output placeholder |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
ab |
`(..., [1]) tensor`
|
Dot product, with shape |
References
- Dot product, Wikipedia.
mdot
Compute the Frobenius inner product of two matrices \(\operatorname{tr}(\mathbf{A}^\mathrm{H}\mathbf{B})\).
Warning
When \(\mathbf{A}\) and \(\mathbf{B}\) are complex, this function returns the complex dot product \((\mathbf{A}, \mathbf{B}) = \operatorname{tr}(\mathbf{A}^\mathrm{H}\mathbf{B})\). This dot product is linear in the second term and antilinear in the first term.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a |
`(..., n, m) tensor`
|
Left matrix, with shape |
required |
b |
`(..., n, m) tensor`
|
Right matrix, with shape |
required |
keepdim |
`bool`
|
Keep reduced dimensions. |
False
|
out |
`(..., [1, 1]) tensor`
|
Output placeholder. |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
dot |
`(..., [1, 1]) tensor`
|
Matrix inner product, with shape |
References
- Frobenius inner product, Wikipedia.
is_orthonormal
Check that a basis is an orthonormal basis.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
basis |
`(F, N, [M]) tensor`
|
A basis of a vector or matrix space, with shape |
required |
return_matrix |
`bool`
|
If True, return the matrix of all pairs of inner products between elements if the basis |
False
|
Returns:
| Name | Type | Description |
|---|---|---|
check |
`bool`
|
True if the basis is orthonormal |
matrix |
`(F, F) tensor`, if `return_matrix is True`
|
Matrix of all pairs of inner products, with shape |