Abstract Class for Special Orthogonal Groups in Vector Representation
Source:R/special-orthogonal.R
SpecialOrthogonalVectors.Rd
Class for the special orthogonal groups
\(\mathrm{SO}(\{2,3\})\) in vector form, i.e. the Lie groups of planar
and 3D rotations. This class is specific to the vector representation of
rotations. For the matrix representation, use the SpecialOrthogonal
class and set n = 2
or n = 3
.
Super classes
rgeomstats::PythonClass
-> rgeomstats::Manifold
-> rgeomstats::LieGroup
-> SpecialOrthogonalVectors
Public fields
n
An integer value specifying the number of rows and columns of the matrices.
epsilon
A numeric value specifying the precision to use for calculations involving potential divison by 0 in rotations.
Methods
Inherited methods
rgeomstats::PythonClass$get_python_class()
rgeomstats::PythonClass$set_python_class()
rgeomstats::Manifold$belongs()
rgeomstats::Manifold$is_tangent()
rgeomstats::Manifold$random_point()
rgeomstats::Manifold$random_tangent_vec()
rgeomstats::Manifold$regularize()
rgeomstats::Manifold$set_metric()
rgeomstats::Manifold$to_tangent()
rgeomstats::LieGroup$add_metric()
rgeomstats::LieGroup$compose()
rgeomstats::LieGroup$exp()
rgeomstats::LieGroup$exp_from_identity()
rgeomstats::LieGroup$exp_not_from_identity()
rgeomstats::LieGroup$get_identity()
rgeomstats::LieGroup$inverse()
rgeomstats::LieGroup$jacobian_translation()
rgeomstats::LieGroup$lie_bracket()
rgeomstats::LieGroup$log()
rgeomstats::LieGroup$log_from_identity()
rgeomstats::LieGroup$log_not_from_identity()
rgeomstats::LieGroup$tangent_translation_map()
Method new()
The SpecialOrthogonalVectors
class constructor.
Usage
SpecialOrthogonalVectors$new(n, epsilon = 0, py_cls = NULL)
Arguments
n
An integer value specifying the number of rows and columns of the matrices.
epsilon
A numeric value specifying the precision to use for calculations involving potential divison by 0 in rotations. Defaults to
0
.py_cls
A Python object of class
SpecialOrthogonalVectors
. Defaults toNULL
in which case it is instantiated on the fly using the other input arguments.
Method projection()
Projects a matrix on \(\mathrm{SO}(2)\) or \(\mathrm{SO}(3)\) using the Frobenius norm.
Arguments
point
A numeric array of shape \([\dots \times n \times n]\) specifying one or more matrices to be projected.
Examples
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$projection(diag(1, 2))
}
Method skew_matrix_from_vector()
Gets the skew-symmetric matrix derived from the vector. In
3D, computes the skew-symmetric matrix, known as the cross-product of a
vector, associated to the vector vec
.
Arguments
vec
A numeric array of shape \([\dots \times \mathrm{dim}]\) specifying one or more vectors from which to compute corresponding skew matrix representations.
Returns
A numeric array of shape \([\dots \times n \times n]\) storing the corresponding skew matrix representations.
Examples
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$skew_matrix_from_vector(array(0))
}
Method vector_from_skew_matrix()
Derives a vector from the skew-symmetric matrix. In 3D, computes the vector defining the cross-product associated to a skew-symmetric matrix.
Returns
A numeric array of shape \([\dots \times \mathrm{dim}]\) storing the corresponding vector representations.
Examples
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$vector_from_skew_matrix(diag(0, 2))
}
Method regularize_tangent_vec_at_identity()
Regularizes a tangent vector at the identity. In 2D, regularizes a tangent vector by getting its norm at the identity to be less than \(\pi\).
Arguments
tangent_vec
A numeric array of shape \([\dots \times 1]\) specifying one or more tangent vectors at base point.
metric
An object of class
RiemannianMetric
specifying the metric to compute the norm of the tangent vector orNULL
. If it is set toNULL
, it defaults to using the Euclidean metric.
Examples
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$regularize_tangent_vec_at_identity(array(0))
}
Method regularize_tangent_vec()
Regularizes a tangent vector at a base point. In 2D, regularizes a tangent vector by getting the norm of its parallel transport to the identity, determined by the metric, to be less than \(\pi\).
Arguments
tangent_vec
A numeric array of shape \([\dots \times 1]\) specifying one or more tangent vectors at corresponding base points.
base_point
A numeric array of shape \([\dots \times 1]\) specifying one or more points on the manifold.
metric
An object of class
RiemannianMetric
specifying the metric to compute the norm of the tangent vector orNULL
. If it is set toNULL
, it defaults to using the Euclidean metric.
Examples
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$regularize_tangent_vec(array(0), array(1))
}
Examples
## ------------------------------------------------
## Method `SpecialOrthogonalVectors$projection`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$projection(diag(1, 2))
}
## ------------------------------------------------
## Method `SpecialOrthogonalVectors$skew_matrix_from_vector`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$skew_matrix_from_vector(array(0))
}
## ------------------------------------------------
## Method `SpecialOrthogonalVectors$vector_from_skew_matrix`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$vector_from_skew_matrix(diag(0, 2))
}
## ------------------------------------------------
## Method `SpecialOrthogonalVectors$regularize_tangent_vec_at_identity`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$regularize_tangent_vec_at_identity(array(0))
}
## ------------------------------------------------
## Method `SpecialOrthogonalVectors$regularize_tangent_vec`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so2 <- SpecialOrthogonal(n = 2, point_type = "vector")
so2$regularize_tangent_vec(array(0), array(1))
}