Class for Lie groups. In this class, point_type ('vector' or
'matrix') will be used to describe the format of the points on the Lie
group. If point_type is 'vector', the format of the inputs is
dimension, where dimension is the dimension of the Lie group. If
point_type is 'matrix', the format of the inputs is c(n, n) where n
is the parameter of \(\mathrm{GL}(n)\) e.g. the amount of rows and
columns of the matrix.
Super classes
rgeomstats::PythonClass -> rgeomstats::Manifold -> LieGroup
Public fields
lie_algebraAn object of class
MatrixLieAlgebraorNULLrepresenting the tangent space at the identity.left_canonical_metricAn object of class
InvariantMetricrepresenting the left invariant metric that corresponds to the Euclidean inner product at the identity.right_canonical_metricAn object of class
InvariantMetricrepresenting the left invariant metric that corresponds to the Euclidean inner product at the identity.metricsA list of objects of class
RiemannianMetric.
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()
Method new()
The LieGroup class constructor.
Usage
LieGroup$new(dim, shape, lie_algebra = NULL, ..., py_cls = NULL)Arguments
dimAn integer value specifying the dimension of the manifold.
shapeAn integer vector specifying the shape of one element of the Lie group.
lie_algebraAn object of class
MatrixLieAlgebraorNULLspecifying the tangent space at the identity....Extra arguments to be passed to parent class constructors. See
Manifoldclass.py_clsA Python object of class
LieGroup. Defaults toNULLin which case it is instantiated on the fly using the other input arguments.
Method exp()
Exponentiates a left-invariant vector field from a base point.
Arguments
tangent_vecA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more tangent vectors at corresponding base points.
base_pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more base points on the manifold. Defaults to identity if
NULL.
Details
The vector input is not an element of the Lie algebra, but of
the tangent space at base_point: if \(g\) denotes base_point,
\(v\) the tangent vector, and \(V = g^{-1} v\) the associated Lie
algebra vector, then $$\exp(v, g) = \mathrm{mul}(g, \exp(V))$$.
Therefore, the Lie exponential is obtained when base_point is NULL,
or the identity.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) storing the group exponential of the input tangent vector(s).
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$exp(rep(0, 3))
}Method exp_from_identity()
Compute the group exponential of tangent vector from the identity.
Arguments
tangent_vecA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more tangent vectors at corresponding base points.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) storing the group exponential of the input tangent vector(s).
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$exp_from_identity(rep(0, 3))
}Method exp_not_from_identity()
Calculate the group exponential at base_point.
Arguments
tangent_vecA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more tangent vectors at corresponding base points.
base_pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more base points on the manifold. Defaults to identity if
NULL.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) storing the group exponential of the input tangent vector(s).
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$exp_not_from_identity(rep(0, 3), rep(0, 3))
}Method log()
Computes a left-invariant vector field bringing base_point
to point.
Arguments
pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more points on the manifold.
base_pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more base points on the manifold. Defaults to identity if
NULL.
Details
The output is a vector of the tangent space at base_point, so
not a Lie algebra element if base_point is not the identity.
Furthermore, denoting point by \(g\) and base_point by \(h\),
the output satisfies $$g = \exp(\log(g, h), h)$$.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) storing the group logarithm of the input point(s).
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$log(rep(0, 3))
}Method log_from_identity()
Computes the group logarithm of point from the identity.
Arguments
pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more points on the manifold.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) storing the group logarithm of the input point(s).
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$log_from_identity(rep(0, 3))
}Method log_not_from_identity()
Computes the group logarithm at base_point.
Arguments
pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more points on the manifold.
base_pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more base points on the manifold. Defaults to identity if
NULL.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) storing the group logarithm of the input point(s).
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$log_not_from_identity(rep(0, 3), rep(0, 3))
}Method get_identity()
Gets the identity of the group.
Returns
A numeric array of shape \(\{ \mathrm{dim}, [n \times n] \}\) storing the identity of the Lie group.
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$get_identity()
}Method lie_bracket()
Computes the lie bracket of two tangent vectors.
Arguments
tangent_vector_aA numeric array of shape \([\dots \times n \times n]\) specifying one or more tangent vectors at corresponding base points.
tangent_vector_bA numeric array of shape \([\dots \times n \times n]\) specifying one or more tangent vectors at corresponding base points.
base_pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more base points on the manifold. Defaults to identity if
NULL.
Details
For matrix Lie groups with tangent vectors \(A\) and \(B\) at the same base point \(P\), this is given by (translate to identity, compute commutator, go back): $$[A,B] = A_P^{-1}B - B_P^{-1}A$$.
Returns
A numeric array of shape \([\dots \times n \times n]\) storing the Lie bracket of the two input tangent vectors.
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$lie_bracket(diag(0, 3), diag(0, 3))
}Method tangent_translation_map()
Computes the push-forward map by the left/right translation.
Arguments
pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more points on the manifold.
left_or_rightA character string specifying whether to compute the map for the left or right translation. Choices are
"left"or"right. Defaults to"left".inverseA boolean specifying whether to inverse the Jacobian matrix. If set to
TRUE, the push forward by the translation by the inverse of the point is returned. Defaults toFALSE.
Details
Computes the push-forward map of the left/right translation by
the point. It corresponds to the tangent map, or differential of the
group multiplication by the point or its inverse. For groups with a
vector representation, it is only implemented at identity, but it can
be used at other points with inverse = TRUE. This method wraps the
Jacobian translation which actually computes the matrix representation
of the map.
Returns
A function computing the tangent map of the left/right
translation by point. It can be applied to tangent vectors.
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$tangent_translation_map(rep(0, 3))
}Method compose()
Performs function composition corresponding to the Lie group.
Arguments
point_aA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more left factors in the product.
point_bA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more right factors in the product.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n
\times n] \}]\) storing the product of point_a and point_b along the
first dimension.
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$compose(rep(0, 3), rep(0, 3))
}Method jacobian_translation()
Computes the Jacobian of left/right translation by a point.
Arguments
pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more points on the manifold.
left_or_rightA character string specifying whether to compute the map for the left or right translation. Choices are
"left"or"right. Defaults to"left".
Returns
A numeric array of shape \([\dots \times \mathrm{dim} \times
\mathrm{dim}]\) storing the Jacobian of the left/right translation by
point.
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$jacobian_translation(rep(0, 3))
}Method inverse()
Computes the inverse law of the Lie group.
Arguments
pointA numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) specifying one or more points to be inverted.
Returns
A numeric array of shape \([\dots \times \{ \mathrm{dim}, [n \times n] \}]\) storing the inverted points.
Examples
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$inverse(rep(0, 3))
}Method add_metric()
Adds a metric to the class $metrics attribute.
Arguments
metricAn object of class
RiemannianMetric.
Examples
## ------------------------------------------------
## Method `LieGroup$exp`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$exp(rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$exp_from_identity`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$exp_from_identity(rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$exp_not_from_identity`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$exp_not_from_identity(rep(0, 3), rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$log`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$log(rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$log_from_identity`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$log_from_identity(rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$log_not_from_identity`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$log_not_from_identity(rep(0, 3), rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$get_identity`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$get_identity()
}
## ------------------------------------------------
## Method `LieGroup$lie_bracket`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$lie_bracket(diag(0, 3), diag(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$tangent_translation_map`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$tangent_translation_map(rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$compose`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$compose(rep(0, 3), rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$jacobian_translation`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$jacobian_translation(rep(0, 3))
}
## ------------------------------------------------
## Method `LieGroup$inverse`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
so3 <- SpecialOrthogonal(n = 3, point_type = "vector")
so3$inverse(rep(0, 3))
}
