Orthogonal and Unitary matrices
Both OrthogonalMatrices
and UnitaryMatrices
are quite similar, as are Rotations
, as well as unitary matrices with determinant equal to one. So these share a {common implementation}(@ref generalunitarymatrices)
Orthogonal Matrices
Manifolds.OrthogonalMatrices
— Type OrthogonalMatrices{n} = GeneralUnitaryMatrices{n,ℝ,AbsoluteDeterminantOneMatrices}
The manifold of (real) orthogonal matrices $\mathrm{O}(n)$.
OrthogonalMatrices(n)
Unitary Matrices
Manifolds.UnitaryMatrices
— Typeconst UnitaryMatrices{n,𝔽} = AbstarctUnitaryMatrices{n,𝔽,AbsoluteDeterminantOneMatrices}
The manifold $U(n,𝔽)$ of $n×n$ complex matrices (when 𝔽=ℂ) or quaternionic matrices (when 𝔽=ℍ) such that
$p^{\mathrm{H}}p = \mathrm{I}_n,$
where $\mathrm{I}_n$ is the $n×n$ identity matrix. Such matrices p
have a property that $\lVert \det(p) \rVert = 1$.
The tangent spaces are given by
\[ T_pU(n) \coloneqq \bigl\{ X \big| pY \text{ where } Y \text{ is skew symmetric, i. e. } Y = -Y^{\mathrm{H}} \bigr\}\]
But note that tangent vectors are represented in the Lie algebra, i.e. just using $Y$ in the representation above.
Constructor
UnitaryMatrices(n, 𝔽::AbstractNumbers=ℂ)
see also OrthogonalMatrices
for the real valued case.
ManifoldsBase.manifold_dimension
— Methodmanifold_dimension(M::UnitaryMatrices{n,ℂ}) where {n}
Return the dimension of the manifold unitary matrices.
\[\dim_{\mathrm{U}(n)} = n^2.\]
ManifoldsBase.manifold_dimension
— Methodmanifold_dimension(M::UnitaryMatrices{n,ℍ})
Return the dimension of the manifold unitary matrices.
\[\dim_{\mathrm{U}(n, ℍ)} = n(2n+1).\]
Common functions
Manifolds.AbsoluteDeterminantOneMatrices
— TypeAbsoluteDeterminantOneMatrices <: AbstractMatrixType
A type to indicate that we require (orthogonal / unitary) matrices with normed determinant, i.e. that the absolute value of the determinant is 1.
Manifolds.AbstractMatrixType
— TypeAbstractMatrixType
A plain type to distinguish different types of matrices, for example DeterminantOneMatrices
and AbsoluteDeterminantOneMatrices
Manifolds.DeterminantOneMatrices
— TypeDeterminantOneMatrices <: AbstractMatrixType
A type to indicate that we require special (orthogonal / unitary) matrices, i.e. of determinant 1.
Manifolds.GeneralUnitaryMatrices
— TypeGeneralUnitaryMatrices{n,𝔽,S<:AbstractMatrixType} <: AbstractDecoratorManifold
A common parametric type for matrices with a unitary property of size $n×n$ over the field $\mathbb F$ which additionally have the AbstractMatrixType
, e.g. are DeterminantOneMatrices
.
Base.exp
— Methodexp(M::Rotations, p, X)
exp(M::OrthogonalMatrices, p, X)
exp(M::UnitaryMatrices, p, X)
Compute the exponential map, that is, since $X$ is represented in the Lie algebra,
exp_p(X) = p\mathrm{e}^X
For different sizes, like $n=2,3,4$ there is specialised implementations
The algorithm used is a more numerically stable form of those proposed in [Gallier2002] and [Andrica2013].
Base.log
— Methodlog(M::Rotations, p, X)
log(M::OrthogonalMatrices, p, X)
log(M::UnitaryMatrices, p, X)
Compute the logarithmic map, that is, since the resulting $X$ is represented in the Lie algebra,
log_p q = \log(p^{\mathrm{H}q)
which is projected onto the skew symmetric matrices for numerical stability.
Base.log
— Methodlog(M::Rotations, p, q)
Compute the logarithmic map on the Rotations
manifold M
which is given by
\[\log_p q = \operatorname{log}(p^{\mathrm{T}}q)\]
where $\operatorname{Log}$ denotes the matrix logarithm. For numerical stability, the result is projected onto the set of skew symmetric matrices.
For antipodal rotations the function returns deterministically one of the tangent vectors that point at q
.
Manifolds.cos_angles_4d_rotation_matrix
— Methodcos_angles_4d_rotation_matrix(R)
4D rotations can be described by two orthogonal planes that are unchanged by the action of the rotation (vectors within a plane rotate only within the plane). The cosines of the two angles $α,β$ of rotation about these planes may be obtained from the distinct real parts of the eigenvalues of the rotation matrix. This function computes these more efficiently by solving the system
\[\begin{aligned} \cos α + \cos β &= \frac{1}{2} \operatorname{tr}(R)\\ \cos α \cos β &= \frac{1}{8} \operatorname{tr}(R)^2 - \frac{1}{16} \operatorname{tr}((R - R^T)^2) - 1. \end{aligned}\]
By convention, the returned values are sorted in decreasing order. See also angles_4d_skew_sym_matrix
.
ManifoldsBase.check_point
— Methodcheck_point(M::Rotations, p; kwargs...)
Check whether p
is a valid point on the UnitaryMatrices
M
, i.e. that $p$ has an determinante of absolute value one, i.e. that $p^{\mathrm{H}}p$
The tolerance for the last test can be set using the kwargs...
.
ManifoldsBase.check_point
— Methodcheck_point(M::UnitaryMatrices, p; kwargs...)
check_point(M::OrthogonalMatrices, p; kwargs...)
check_point(M::GeneralUnitaryMatrices{n,𝔽}, p; kwargs...)
Check whether p
is a valid point on the UnitaryMatrices
or [OrthogonalMatrices
] M
, i.e. that $p$ has an determinante of absolute value one
The tolerance for the last test can be set using the kwargs...
.
ManifoldsBase.check_vector
— Methodcheck_vector(M::UnitaryMatrices{n}, p, X; kwargs... )
check_vector(M::OrthogonalMatrices{n}, p, X; kwargs... )
check_vector(M::Rotations{n}, p, X; kwargs... )
check_vector(M::GeneralUnitaryMatrices{n,𝔽}, p, X; kwargs... )
Check whether X
is a tangent vector to p
on the UnitaryMatrices
space M
, i.e. after check_point
(M,p)
, X
has to be skew symmetric (Hermitian) and orthogonal to p
.
The tolerance for the last test can be set using the kwargs...
.
ManifoldsBase.embed
— Methodembed(M::GeneralUnitaryMatrices{n,𝔽}, p, X)
Embed the tangent vector X
at point p
in M
from its Lie algebra representation (set of skew matrices) into the Riemannian submanifold representation
The formula reads
\[X_{\text{embedded}} = p * X\]
ManifoldsBase.get_coordinates
— Methodget_coordinates(M::Rotations, p, X)
get_coordinates(M::OrthogonalMatrices, p, X)
get_coordinates(M::UnitaryMatrices, p, X)
Extract the unique tangent vector components $X^i$ at point p
on Rotations
$\mathrm{SO}(n)$ from the matrix representation X
of the tangent vector.
The basis on the Lie algebra $𝔰𝔬(n)$ is chosen such that for $\mathrm{SO}(2)$, $X^1 = θ = X_{21}$ is the angle of rotation, and for $\mathrm{SO}(3)$, $(X^1, X^2, X^3) = (X_{32}, X_{13}, X_{21}) = θ u$ is the angular velocity and axis-angle representation, where $u$ is the unit vector along the axis of rotation.
For $\mathrm{SO}(n)$ where $n ≥ 4$, the additional elements of $X^i$ are $X^{j (j - 3)/2 + k + 1} = X_{jk}$, for $j ∈ [4,n], k ∈ [1,j)$.
ManifoldsBase.get_embedding
— Methodget_embedding(M::OrthogonalMatrices{n})
get_embedding(M::Rotations{n})
get_embedding(M::UnitaryMatrices{n})
Return the embedding, i.e. The $\mathbb F^{n×n}$, where $\mathbb F = \mathbb R$ for the first two and $\mathbb F = \mathbb C$ for the unitary matrices.
ManifoldsBase.get_vector
— Methodget_vector(M::OrthogonalMatrices, p, Xⁱ, B::DefaultOrthogonalBasis)
get_vector(M::Rotations, p, Xⁱ, B::DefaultOrthogonalBasis)
Convert the unique tangent vector components Xⁱ
at point p
on Rotations
or OrthogonalMatrices
to the matrix representation $X$ of the tangent vector. See get_coordinates
for the conventions used.
ManifoldsBase.injectivity_radius
— Methodinjectivity_radius(G::GeneraliUnitaryMatrices)
Return the injectivity radius for general unitary matrix manifolds, which is[1]
\[ \operatorname{inj}_{\mathrm{U}(n)} = π.\]
ManifoldsBase.injectivity_radius
— Methodinjectivity_radius(G::GeneralUnitaryMatrices{n,ℂ,DeterminantOneMatrices})
Return the injectivity radius for general complex unitary matrix manifolds, where the determinant is $+1$, which is[1]
\[ \operatorname{inj}_{\mathrm{SU}(n)} = π \sqrt{2}.\]
ManifoldsBase.injectivity_radius
— Methodinjectivity_radius(G::SpecialOrthogonal)
injectivity_radius(G::Orthogonal)
injectivity_radius(M::Rotations)
injectivity_radius(M::Rotations, ::ExponentialRetraction)
Return the radius of injectivity on the Rotations
manifold M
, which is $π\sqrt{2}$. [1]
ManifoldsBase.is_flat
— Methodis_flat(M::GeneralUnitaryMatrices)
Return true if GeneralUnitaryMatrices
M
is SO(2) or U(1) and false otherwise.
ManifoldsBase.manifold_dimension
— Methodmanifold_dimension(M::GeneralUnitaryMatrices{n,ℂ,DeterminantOneMatrices})
Return the dimension of the manifold of special unitary matrices.
\[\dim_{\mathrm{SU}(n)} = n^2-1.\]
ManifoldsBase.manifold_dimension
— Methodmanifold_dimension(M::Rotations)
manifold_dimension(M::OrthogonalMatrices)
Return the dimension of the manifold orthogonal matrices and of the manifold of rotations
\[\dim_{\mathrm{O}(n)} = \dim_{\mathrm{SO}(n)} = \frac{n(n-1)}{2}.\]
ManifoldsBase.project
— Method project(M::OrthogonalMatrices{n}, p, X)
project(M::Rotations{n}, p, X)
project(M::UnitaryMatrices{n}, p, X)
Orthogonally project the tangent vector $X ∈ 𝔽^{n × n}$, $\mathbb F ∈ \{\mathbb R, \mathbb C\}$ to the tangent space of M
at p
, and change the representer to use the corresponding Lie algebra, i.e. we compute
\[ \operatorname{proj}_p(X) = \frac{p^{\mathrm{H}} X - (p^{\mathrm{H}} X)^{\mathrm{H}}}{2},\]
ManifoldsBase.project
— Method project(G::UnitaryMatrices{n}, p)
project(G::OrthogonalMatrices{n}, p)
Project the point $p ∈ 𝔽^{n × n}$ to the nearest point in $\mathrm{U}(n,𝔽)=$Unitary(n,𝔽)
under the Frobenius norm. If $p = U S V^\mathrm{H}$ is the singular value decomposition of $p$, then the projection is
\[ \operatorname{proj}_{\mathrm{U}(n,𝔽)} \colon p ↦ U V^\mathrm{H}.\]
ManifoldsBase.retract
— Methodretract(M::Rotations, p, X, ::PolarRetraction)
retract(M::OrthogonalMatrices, p, X, ::PolarRetraction)
Compute the SVD-based retraction on the Rotations
and OrthogonalMatrices
M
from p
in direction X
(as an element of the Lie group) and is a second-order approximation of the exponential map. Let
\[USV = p + pX\]
be the singular value decomposition, then the formula reads
\[\operatorname{retr}_p X = UV^\mathrm{T}.\]
ManifoldsBase.retract
— Methodretract(M::Rotations, p, X, ::QRRetraction)
retract(M::OrthogonalMatrices, p. X, ::QRRetraction)
Compute the QR-based retraction on the Rotations
and OrthogonalMatrices
M
from p
in direction X
(as an element of the Lie group), which is a first-order approximation of the exponential map.
This is also the default retraction on these manifolds.
ManifoldsBase.riemann_tensor
— Methodriemann_tensor(::GeneralUnitaryMatrices, p, X, Y, Z)
Compute the value of Riemann tensor on the GeneralUnitaryMatrices
manifold. The formula reads[Rentmeesters2011] $R(X,Y)Z=\frac{1}{4}[Z, [X, Y]]$.
Statistics.mean
— Methodmean(
M::Rotations,
x::AbstractVector,
[w::AbstractWeights,]
method = GeodesicInterpolationWithinRadius(π/2/√2);
kwargs...,
)
Compute the Riemannian mean
of x
using GeodesicInterpolationWithinRadius
.
Footnotes and References
- Gallier2002
Gallier J.; Xu D.; Computing exponentials of skew-symmetric matrices and logarithms of orthogonal matrices. International Journal of Robotics and Automation (2002), 17(4), pp. 1-11. pdf.
- Andrica2013
Andrica D.; Rohan R.-A.; Computing the Rodrigues coefficients of the exponential map of the Lie groups of matrices. Balkan Journal of Geometry and Its Applications (2013), 18(2), pp. 1-2. pdf.
- 1
For a derivation of the injectivity radius, see sethaxen.com/blog/2023/02/the-injectivity-radii-of-the-unitary-groups/.
- Rentmeesters2011
Q. Rentmeesters, “A gradient method for geodesic data fitting on some symmetric Riemannian manifolds,” in 2011 50th IEEE Conference on Decision and Control and European Control Conference, Dec. 2011, pp. 7141–7146. doi: 10.1109/CDC.2011.6161280.