Rotations
The manifold $\mathrm{SO}(n)$ of orthogonal matrices with determinant $+1$ in $ℝ^{n × n}$, i.e.
The Lie group $\mathrm{SO}(n)$ is a subgroup of the orthogonal group $\mathrm{O}(n)$ and also known as the special orthogonal group or the set of rotations group. See also SpecialOrthogonal
, which is this manifold equipped with the group operation.
Tangent vectors are represented by elements of the corresponding Lie algebra, which is the tangent space at the identity element. This convention allows for more efficient operations on tangent vectors. Tangent spaces at different points are different vector spaces.
Let $L_R: \mathrm{SO}(n) → \mathrm{SO}(n)$ where $R ∈ \mathrm{SO}(n)$ be the left-multiplication by $R$, that is $L_R(S) = RS$. The tangent space at rotation $R$, $T_R \mathrm{SO}(n)$, is related to the tangent space at the identity rotation $I_n$ by the differential of $L_R$ at identity, $(\mathrm{d}L_R)_{I_n} : T_{I_n} \mathrm{SO}(n) → T_R \mathrm{SO}(n)$. For a tangent vector at the identity rotation $X ∈ T_{I_n} \mathrm{SO}(n)$ the matrix representation of the corresponding tangent vector $Y$ at a rotation $R$ can be obtained by matrix multiplication: $Y = RX ∈ T_R \mathrm{SO}(n)$. You can compare the functions log
and exp
to see how it works in practice.
Manifolds.NormalRotationDistribution
— TypeNormalRotationDistribution(M::Rotations, d::Distribution, x::TResult)
Distribution that returns a random point on the manifold Rotations
M
. Random point is generated using base distribution d
and the type of the result is adjusted to TResult
.
See normal_rotation_distribution
for details.
Manifolds.Rotations
— TypeRotations{N} <: Manifold{ℝ}
The special orthogonal manifold $\mathrm{SO}(n)$ represented by $n × n$ real-valued orthogonal matrices with determinant $+1$ is the manifold of Rotations
, since these matrices represent all rotations of points in $ℝ^n$.
Constructor
Rotations(n)
Generate the $\mathrm{SO}(n) \subset ℝ^{n × n}$
Base.exp
— Methodexp(M::Rotations, p, X)
Compute the exponential map on the Rotations
from p
into direction X
, i.e.
where $\operatorname{Exp}(X)$ denotes the matrix exponential of $X$.
exp(M::Rotations{4}, p, X)
Compute the exponential map of tangent vector X
at point p
from $\mathrm{SO}(4)$ manifold M
.
The algorithm used is a more numerically stable form of those proposed in [Gallier2002] and [Andrica2013].
Base.log
— Methodlog(M::Rotations, p, q)
Compute the logarithmic map on the Rotations
manifold M
$=\mathrm{SO}(n)$, which is given by
where $\operatorname{Log}$ denotes the matrix logarithm.
For antipodal rotations the function returns deterministically one of the tangent vectors that point at q
.
LinearAlgebra.norm
— Methodnorm(M::Rotations, p, X)
Compute the norm of a tangent vector X
from the tangent space at p
on the Rotations
M
. The formula reads
i.e. the Frobenius norm of X
, where tangent vectors are represented by elements from the Lie algebra.
Manifolds.angles_4d_skew_sym_matrix
— Methodangles_4d_skew_sym_matrix(A)
The Lie algebra of Rotations(4)
in $ℝ^{4 × 4}$, $𝔰𝔬(4)$, consists of $4 × 4$ skew-symmetric matrices. The unique imaginary components of their eigenvalues are the angles of the two plane rotations. This function computes these more efficiently than eigvals
.
By convention, the returned values are sorted in decreasing order (corresponding to the same ordering of angles as cos_angles_4d_rotation_matrix
).
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
By convention, the returned values are sorted in increasing order. See angles_4d_skew_sym_matrix
.
Manifolds.normal_rotation_distribution
— Methodnormal_rotation_distribution(M::Rotations, p, σ::Real)
Return a random point on the manifold Rotations
M
by generating a (Gaussian) random orthogonal matrix with determinant $+1$. Let
be the QR decomposition of a random matrix $A$, then the formula reads
where $D$ is a diagonal matrix with the signs of the diagonal entries of $R$, i.e.
It can happen that the matrix gets -1 as a determinant. In this case, the first and second columns are swapped.
The argument p
is used to determine the type of returned points.
Manifolds.normal_tvector_distribution
— Methodnormal_tvector_distribution(M::Rotations, p, σ)
Normal distribution in ambient space with standard deviation σ
projected to tangent space at p
.
ManifoldsBase.check_manifold_point
— Methodcheck_manifold_point(M, p; kwargs...)
Check whether p
is a valid point on the Rotations
M
, i.e. is an array of size manifold_dimension
(M)
and represents a valid rotation. The tolerance for the last test can be set using the kwargs...
.
ManifoldsBase.check_tangent_vector
— Methodcheck_tangent_vector(M, p, X; check_base_point = true, kwargs... )
Check whether X
is a tangent vector to p
on the Rotations
space M
, i.e. after check_manifold_point
(M,p)
, X
has to be of same dimension and orthogonal to p
. The optional parameter check_base_point
indicates, whether to call check_manifold_point
for p
. The tolerance for the last test can be set using the kwargs...
.
ManifoldsBase.get_coordinates
— Methodget_coordinates(M::Rotations, 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_vector
— Methodget_vector(M::Rotations, p, Xⁱ, B:: DefaultOrthogonalBasis)
Convert the unique tangent vector components Xⁱ
at point p
on Rotations
group $\mathrm{SO}(n)$ to the matrix representation $X$ of the tangent vector. See get_coordinates
for the conventions used.
ManifoldsBase.injectivity_radius
— Methodinjectivity_radius(M::Rotations)
injectivity_radius(M::Rotations, p)
Return the injectivity radius on the Rotations
M
, which is globally
injectivity_radius(M::Rotations, p, ::PolarRetraction)
Return the radius of injectivity for the PolarRetraction
on the Rotations
M
which is $\frac{π}{\sqrt{2}}$.
ManifoldsBase.inner
— Methodinner(M::Rotations, p, X, Y)
Compute the inner product of the two tangent vectors X
, Y
from the tangent plane at p
on the special orthogonal space M=
$\mathrm{SO}(n)$ using the restriction of the metric from the embedding, i.e.
Tangent vectors are represented by matrices.
ManifoldsBase.inverse_retract
— Methodinverse_retract(M, p, q, ::PolarInverseRetraction)
Compute a vector from the tangent space $T_p\mathrm{SO}(n)$ of the point p
on the Rotations
manifold M
with which the point q
can be reached by the PolarRetraction
from the point p
after time 1.
The formula reads
where $s$ is the solution to the Sylvester equation
ManifoldsBase.inverse_retract
— Methodinverse_retract(M::Rotations, p, q, ::QRInverseRetraction)
Compute a vector from the tangent space $T_p\mathrm{SO}(n)$ of the point p
on the Rotations
manifold M
with which the point q
can be reached by the QRRetraction
from the point q
after time 1.
ManifoldsBase.manifold_dimension
— Methodmanifold_dimension(M::Rotations)
Return the dimension of the manifold $\mathrm{SO}(n)$, i.e.
ManifoldsBase.project
— Methodproject(M::Rotations, p, X)
Project the matrix X
onto the tangent space by making X
skew symmetric,
where tangent vectors are represented by elements from the Lie group
ManifoldsBase.project
— Methodproject(M::Rotations, p; check_det = true)
Project p
to the nearest point on manifold M
.
Given the singular value decomposition $p = U Σ V^\mathrm{T}$, with the singular values sorted in descending order, the projection is
The diagonal matrix ensures that the determinant of the result is $+1$. If p
is expected to be almost special orthogonal, then you may avoid this check with check_det = false
.
ManifoldsBase.representation_size
— Methodrepresentation_size(M::Rotations)
Return the size()
of a point on the Rotations
M
, i.e. for the $\mathrm{SO}(n)$ it's (n,n)
.
ManifoldsBase.retract
— Methodretract(M::Rotations, p, X, ::PolarRetraction)
Compute the SVD-based retraction on the Rotations
M
from p
in direction X
(as an element of the Lie group) and is a second-order approximation of the exponential map. Let
be the singular value decomposition, then the formula reads
ManifoldsBase.retract
— Methodretract(M, p, X, ::QRRetraction)
Compute the QR-based retraction on the Rotations
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 the Rotations
ManifoldsBase.zero_tangent_vector
— Methodzero_tangent_vector(M::Rotations, p)
Return the zero tangent vector from the tangent space art p
on the Rotations
as an element of the Lie group, i.e. the zero matrix.
Statistics.mean
— Methodmean(
M::Rotations,
x::AbstractVector,
[w::AbstractWeights,]
method = GeodesicInterpolationWithinRadius(π/2/√2);
kwargs...,
)
Compute the Riemannian mean
of x
using GeodesicInterpolationWithinRadius
.
Literature
- 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.