An R6::R6Class object implementing the base
RiemannianMetric
class. This is an abstract class for Riemannian and
pseudo-Riemannian metrics which are the associated Levi-Civita connection
on the tangent bundle.
Super classes
rgeomstats::PythonClass
-> rgeomstats::Connection
-> RiemannianMetric
Methods
Inherited methods
rgeomstats::PythonClass$get_python_class()
rgeomstats::PythonClass$set_python_class()
rgeomstats::Connection$christoffels()
rgeomstats::Connection$curvature()
rgeomstats::Connection$curvature_derivative()
rgeomstats::Connection$directional_curvature()
rgeomstats::Connection$directional_curvature_derivative()
rgeomstats::Connection$exp()
rgeomstats::Connection$geodesic()
rgeomstats::Connection$geodesic_equation()
rgeomstats::Connection$injectivity_radius()
rgeomstats::Connection$ladder_parallel_transport()
rgeomstats::Connection$log()
rgeomstats::Connection$parallel_transport()
Method new()
The RiemannianMetric
class constructor.
Usage
RiemannianMetric$new(
dim,
shape = NULL,
signature = NULL,
default_coords_type = "intrinsic",
py_cls = NULL
)
Arguments
dim
An integer value specifying the dimension of the manifold.
shape
An integer vector specifying the shape of one element of the manifold. Defaults to
NULL
.signature
An integer vector specifying the signature of the metric. Defaults to
c(dim, 0L)
.default_coords_type
A string specifying the coordinate type. Choices are
extrensic
orintrinsic
. Defaults tointrinsic
.py_cls
A Python object of class
RiemannianMetric
. Defaults toNULL
in which case it is instantiated on the fly using the other input arguments.
Method metric_matrix()
Metric matrix at the tangent space at a base point.
Arguments
base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.
Examples
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$metric_matrix()
}
Method cometric_matrix()
Inner co-product matrix at the cotangent space at a base point. This represents the cometric matrix, i.e. the inverse of the metric matrix.
Arguments
base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.
Examples
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$cometric_matrix()
}
Method inner_product_derivative_matrix()
Compute derivative of the inner prod matrix at base point.
Arguments
base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.
Returns
A numeric array of shape dim x dim
storing the derivative of
the inverse of the inner-product matrix.
Examples
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$inner_product_derivative_matrix()
}
Method inner_product()
Inner product between two tangent vectors at a base point.
Arguments
tangent_vec_a
A numeric array of shape
dim
specifying a tangent vector at base point.tangent_vec_b
A numeric array of shape
dim
specifying a tangent vector at base point.base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.
Returns
A scalar value representing the inner product between the two input tangent vectors at the input base point.
Examples
if (reticulate::py_module_available("geomstats")) {
mt <- SPDMetricLogEuclidean$new(n = 3)
mt$inner_product(diag(0, 3), diag(1, 3), base_point = diag(1, 3))
}
Method inner_coproduct()
Computes inner coproduct between two cotangent vectors at base point. This is the inner product associated to the cometric matrix.
Arguments
cotangent_vec_a
A numeric array of shape
dim
specifying a cotangent vector at base point.cotangent_vec_b
A numeric array of shape
dim
specifying a cotangent vector at base point.base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.
Returns
A scalar value representing the inner coproduct between the two input cotangent vectors at the input base point.
Examples
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$inner_coproduct(diag(0, 3), diag(1, 3), base_point = diag(1, 3))
}
Method hamiltonian()
Computes the Hamiltonian energy associated to the cometric. The Hamiltonian at state \((q, p)\) is defined by $$H(q, p) = \frac{1}{2} \langle p, p \rangle_q,$$ where \(\langle \cdot, \cdot \rangle_q\) is the cometric at \(q\).
Arguments
state
A list with two components: (i) a numeric array of shape
dim
specifying the position which is a point on the manifold and (ii) a numeric array of shapedim
specifying the momentum which is a cotangent vector.
Examples
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$hamiltonian()
}
Method squared_norm()
Computes the square of the norm of a vector. Squared norm of a vector associated to the inner product at the tangent space at a base point.
Method norm()
Computes the norm of a vector associated to the inner product at the tangent space at a base point.
Method normalize()
Normalizes a tangent vector at a given point.
Arguments
vector
A numeric array of shape
dim
specifying a vector.base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.
Examples
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$normalize(diag(2, 3), diag(1, 3))
}
Method random_unit_tangent_vec()
Generates a random unit tangent vector at a given point.
Arguments
base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.n_vectors
An integer value specifying the number of vectors to be generated at
base_point
. For vectorization purposes,n_vectors
can be greater than 1 iffbase_point
corresponds to a single point. Defaults to1L
.
Returns
A numeric array of shape c(n_vectors, dim)
storing random unit
tangent vectors at base_point
.
Examples
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$random_unit_tangent_vec(diag(1, 3))
}
Method squared_dist()
Squared geodesic distance between two points.
Arguments
point_a
A numeric array of shape
dim
on the manifold.point_b
A numeric array of shape
dim
on the manifold....
Extra parameters to be passed to the
$log()
method of the parentConnection
class.
Method dist()
Geodesic distance between two points.
Arguments
point_a
A numeric array of shape
dim
on the manifold.point_b
A numeric array of shape
dim
on the manifold....
Extra parameters to be passed to the
$log()
method of the parentConnection
class.
Examples
if (reticulate::py_module_available("geomstats")) {
mt <- SPDMetricLogEuclidean$new(n = 3)
mt$dist(diag(1, 3), diag(1, 3))
}
Method dist_broadcast()
Computes the geodesic distance between points.
Arguments
points_a
A numeric array of shape
c(n_samples_a, dim)
specifying a set of points on the manifold.points_b
A numeric array of shape
c(n_samples_b, dim)
specifying a set of points on the manifold.
Details
If n_samples_a == n_samples_b
then dist
is the element-wise
distance result of a point in points_a
with the point from points_b
of the same index. If n_samples_a != n_samples_b
then dist
is the
result of applying geodesic distance for each point from points_a
to
all points from points_b
.
Returns
A numeric array of shape c(n_samples_a, dim)
if n_samples_a == n_samples_b
or of shape c(n_samples_a, n_samples_b, dim)
if
n_samples_a != n_samples_b
storing the geodesic distance between
points in set A and points in set B.
Examples
if (reticulate::py_module_available("geomstats")) {
mt <- SPDMetricLogEuclidean$new(n = 3)
mt$dist(diag(1, 3), diag(1, 3))
}
Method dist_pairwise()
Computes the pairwise distance between points.
Arguments
points
A numeric array of shape
c(n_samples, dim)
specifying a set of points on the manifold.n_jobs
An integer value specifying the number of cores for parallel computation. Defaults to
1L
....
Extra parameters to be passed tothe
joblib.Parallel
Python class. See joblib documentation for details.
Method diameter()
Computes the diameter of set of points on a manifold.
Method closest_neighbor_index()
Finds the closest neighbor to a point among a set of neighbors.
Method normal_basis()
Normalizes the basis with respect to the metric. This corresponds to a renormalization of each basis vector.
Method sectional_curvature()
Computes the sectional curvature.
Arguments
tangent_vec_a
A numeric array of shape
c(n, n)
specifying a tangent vector atbase_point
.tangent_vec_b
A numeric array of shape
c(n, n)
specifying a tangent vector atbase_point
.base_point
A numeric array of shape
dim
specifying a point on the manifold. Defaults toNULL
.
Details
For two orthonormal tangent vectors \(x\) and \(y\) at a base point, the sectional curvature is defined by $$\langle R(x, y)x, y \rangle = \langle R_x(y), y \rangle.$$ For non-orthonormal vectors, it is $$\langle R(x, y)x, y \rangle / \\|x \wedge y\\|^2.$$ See also https://en.wikipedia.org/wiki/Sectional_curvature.
Examples
## ------------------------------------------------
## Method `RiemannianMetric$metric_matrix`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$metric_matrix()
}
## ------------------------------------------------
## Method `RiemannianMetric$cometric_matrix`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$cometric_matrix()
}
## ------------------------------------------------
## Method `RiemannianMetric$inner_product_derivative_matrix`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$inner_product_derivative_matrix()
}
## ------------------------------------------------
## Method `RiemannianMetric$inner_product`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
mt <- SPDMetricLogEuclidean$new(n = 3)
mt$inner_product(diag(0, 3), diag(1, 3), base_point = diag(1, 3))
}
## ------------------------------------------------
## Method `RiemannianMetric$inner_coproduct`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$inner_coproduct(diag(0, 3), diag(1, 3), base_point = diag(1, 3))
}
## ------------------------------------------------
## Method `RiemannianMetric$hamiltonian`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$hamiltonian()
}
## ------------------------------------------------
## Method `RiemannianMetric$normalize`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$normalize(diag(2, 3), diag(1, 3))
}
## ------------------------------------------------
## Method `RiemannianMetric$random_unit_tangent_vec`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
# mt <- SPDMetricLogEuclidean$new(n = 3)
# mt$random_unit_tangent_vec(diag(1, 3))
}
## ------------------------------------------------
## Method `RiemannianMetric$dist`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
mt <- SPDMetricLogEuclidean$new(n = 3)
mt$dist(diag(1, 3), diag(1, 3))
}
## ------------------------------------------------
## Method `RiemannianMetric$dist_broadcast`
## ------------------------------------------------
if (reticulate::py_module_available("geomstats")) {
mt <- SPDMetricLogEuclidean$new(n = 3)
mt$dist(diag(1, 3), diag(1, 3))
}