Metric manifold
A Riemannian manifold always consists of a topological manifold together with a smoothly varying metric $g$.
However, often there is an implicitly assumed (default) metric, like the usual inner product on Euclidean
space. This decorator takes this into account. It is not necessary to use this decorator if you implement just one (or the first) metric. If you later introduce a second, the old (first) metric can be used with the (non MetricManifold
) Manifold
, i.e. without an explicitly stated metric.
This manifold decorator serves two purposes:
- to implement different metrics (e.g. in closed form) for one
Manifold
- to provide a way to compute geodesics on manifolds, where this
Metric
does not yield closed formula.
Let's first look at the provided types.
Types
Manifolds.Metric
β TypeMetric
Abstract type for the pseudo-Riemannian metric tensor $g$, a family of smoothly varying inner products on the tangent space. See inner
.
Manifolds.MetricManifold
β TypeMetricManifold{π½,M<:Manifold{π½},G<:Metric} <: AbstractDecoratorManifold{π½}
Equip a Manifold
explicitly with a Metric
G
.
For a Metric Manifold, by default, assumes, that you implement the linear form from local_metric
in order to evaluate the exponential map.
If the corresponding Metric
G
yields closed form formulae for e.g. the exponential map and this is implemented directly (without solving the ode), you can of course still implement that directly.
Constructor
MetricManifold(M, G)
Manifolds.RiemannianMetric
β TypeRiemannianMetric <: Metric
Abstract type for Riemannian metrics, a family of positive definite inner products. The positive definite property means that for $X β T_p \mathcal M$, the inner product $g(X, X) > 0$ whenever $X$ is not the zero vector.
Implement Different Metrics on the same Manifold
In order to distinguish different metrics on one manifold, one can introduce two Metric
s and use this type to dispatch on the metric, see SymmetricPositiveDefinite
. To avoid overhead, one Metric
can then be marked as being the default, i.e. the one that is used, when no MetricManifold
decorator is present. This avoids reimplementation of the first existing metric, access to the metric-dependent functions that were implemented using the undecorated manifold, as well as the transparent fallback of the corresponding MetricManifold
with default metric to the undecorated implementations. This does not cause any runtime overhead. Introducing a default Metric
serves a better readability of the code when working with different metrics.
Implementation of Metrics
For the case that a local_metric
is implemented as a bilinear form that is positive definite, the following further functions are provided, unless the corresponding Metric
is marked as default β then the fallbacks mentioned in the last section are used for e.g. the exp!
onential map.
Base.exp
β Methodexp(N::MetricManifold{M,G}, p, X)
Copute the exponential map on the Manifold
M
equipped with the Metric
G
.
If the metric was declared the default metric using is_default_metric
, this method falls back to exp(M,p,X)
.
Otherwise it numerically integrates the underlying ODE, see solve_exp_ode
. Currently, the numerical integration is only accurate when using a single coordinate chart that covers the entire manifold. This excludes coordinates in an embedded space.
Base.log
β Methodlog(N::MetricManifold{M,G}, p, q)
Copute the logarithmic map on the Manifold
M
equipped with the Metric
G
.
If the metric was declared the default metric using is_default_metric
, this method falls back to log(M,p,q)
. Otherwise, you have to provide an implementation for the non-default Metric
G
metric within its MetricManifold
{M,G}
.
Manifolds.christoffel_symbols_first
β Methodchristoffel_symbols_first(
M::MetricManifold,
p;
backend::AbstractDiffBackend = diff_backend(),
)
Compute the Christoffel symbols of the first kind in local coordinates. The Christoffel symbols are (in Einstein summation convention)
where $g_{ij,k}=\frac{β}{β p^k} g_{ij}$ is the coordinate derivative of the local representation of the metric tensor. The dimensions of the resulting multi-dimensional array are ordered $(i,j,k)$.
Manifolds.christoffel_symbols_second
β Methodchristoffel_symbols_second(
M::MetricManifold,
x;
backend::AbstractDiffBackend = diff_backend(),
)
Compute the Christoffel symbols of the second kind in local coordinates. The Christoffel symbols are (in Einstein summation convention)
where $Ξ_{ijk}$ are the Christoffel symbols of the first kind, and $g^{kl}$ is the inverse of the local representation of the metric tensor. The dimensions of the resulting multi-dimensional array are ordered $(l,i,j)$.
Manifolds.christoffel_symbols_second_jacobian
β Methodchristoffel_symbols_second_jacobian(
M::MetricManifold,
p;
backend::AbstractDiffBackend = diff_backend(),
)
Get partial derivatives of the Christoffel symbols of the second kind for manifold M
at p
with respect to the coordinates of p
, $\frac{β}{β p^l} Ξ^{k}_{ij} = Ξ^{k}_{ij,l}.$ The dimensions of the resulting multi-dimensional array are ordered $(i,j,k,l)$.
Manifolds.det_local_metric
β Methoddet_local_metric(M::MetricManifold, p)
Return the determinant of local matrix representation of the metric tensor $g$.
Manifolds.einstein_tensor
β Methodeinstein_tensor(M::MetricManifold, p; backend::AbstractDiffBackend = diff_backend())
Compute the Einstein tensor of the manifold M
at the point p
.
Manifolds.flat
β Methodflat(N::MetricManifold{M,G}, p, X::FVector{TangentSpaceType})
Compute the musical isomorphism to transform the tangent vector X
from the Manifold
M
equipped with Metric
G
to a cotangent by computing
where $G_p$ is the local matrix representation of G
, see local_metric
Manifolds.gaussian_curvature
β Methodgaussian_curvature(M::MetricManifold, x; backend::AbstractDiffBackend = diff_backend())
Compute the Gaussian curvature of the manifold M
at the point x
.
Manifolds.inverse_local_metric
β Methodinverse_local_metric(M::MetricManifold, p)
Return the local matrix representation of the inverse metric (cometric) tensor, usually written $g^{ij}$.
Manifolds.is_default_metric
β Methodis_default_metric(M,G)
Indicate whether the Metric
G
is the default metric for the Manifold
M
. This means that any occurence of MetricManifold
(M,G) where typeof(is_default_metric(M,G)) = true
falls back to just be called with M
such that the Manifold
M
implicitly has this metric, for example if this was the first one implemented or is the one most commonly assumed to be used.
Manifolds.is_default_metric
β Methodis_default_metric(MM::MetricManifold)
Indicate whether the Metric
MM.G
is the default metric for the Manifold
MM.manifold,
within the MetricManifold
MM
. This means that any occurence of MetricManifold
(MM.manifold, MM.G)
where is_default_metric(MM.manifold, MM.G)) = true
falls back to just be called with MM.manifold,
such that the Manifold
MM.manifold
implicitly has the metric MM.G
, for example if this was the first one implemented or is the one most commonly assumed to be used.
Manifolds.local_metric
β Methodlocal_metric(M::MetricManifold, p)
Return the local matrix representation at the point p
of the metric tensor $g$ on the Manifold
M
, usually written $g_{ij}$. The matrix has the property that $g(X, Y)=X^\mathrm{T} [g_{ij}] Y = g_{ij} X^i Y^j$, where the latter expression uses Einstein summation convention.
Manifolds.local_metric_jacobian
β Methodlocal_metric_jacobian(
M::MetricManifold,
p;
backend::AbstractDiffBackend = diff_backend(),
)
Get partial derivatives of the local metric of M
at p
with respect to the coordinates of p
, $\frac{β}{β p^k} g_{ij} = g_{ij,k}$. The dimensions of the resulting multi-dimensional array are ordered $(i,j,k)$.
Manifolds.log_local_metric_density
β Methodlog_local_metric_density(M::MetricManifold, p)
Return the natural logarithm of the metric density $Ο$ of M
at p
, which is given by $Ο = \log \sqrt{|\det [g_{ij}]|}$.
Manifolds.metric
β Methodmetric(M::MetricManifold)
Get the metric $g$ of the manifold M
.
Manifolds.ricci_curvature
β Methodricci_curvature(M::MetricManifold, p; backend::AbstractDiffBackend = diff_backend())
Compute the Ricci scalar curvature of the manifold M
at the point p
.
Manifolds.ricci_tensor
β Methodricci_tensor(M::MetricManifold, p; backend::AbstractDiffBackend = diff_backend())
Compute the Ricci tensor, also known as the Ricci curvature tensor, of the manifold M
at the point p
.
Manifolds.riemann_tensor
β Methodriemann_tensor(M::MetricManifold, p; backend::AbstractDiffBackend = diff_backend())
Compute the Riemann tensor $R^l_{ijk}$, also known as the Riemann curvature tensor, at the point p
. The dimensions of the resulting multi-dimensional array are ordered $(l,i,j,k)$.
Manifolds.sharp
β Methodsharp(N::MetricManifold{M,G}, p, ΞΎ::FVector{CotangentSpaceType})
Compute the musical isomorphism to transform the cotangent vector ΞΎ
from the Manifold
M
equipped with Metric
G
to a tangent by computing
where $G_p$ is the local matrix representation of G
, i.e. one employs inverse_local_metric
here to obtain $G_p^{-1}$.
Manifolds.solve_exp_ode
β Methodsolve_exp_ode(
M::MetricManifold,
p,
X,
tspan;
backend::AbstractDiffBackend = diff_backend(),
solver = AutoVern9(Rodas5()),
kwargs...,
)
Approximate the exponential map on the manifold over the provided timespan assuming the Levi-Civita connection by solving the ordinary differential equation
where $Ξ^k_{ij}$ are the Christoffel symbols of the second kind, and the Einstein summation convention is assumed. The arguments tspan
and solver
follow the OrdinaryDiffEq
conventions. kwargs...
specify keyword arguments that will be passed to OrdinaryDiffEq.solve
.
Currently, the numerical integration is only accurate when using a single coordinate chart that covers the entire manifold. This excludes coordinates in an embedded space.
This function only works for Julia 1.1 or greater, when OrdinaryDiffEq.jl is loaded with
using OrdinaryDiffEq
ManifoldsBase.inner
β Methodinner(N::MetricManifold{M,G}, p, X, Y)
Compute the inner product of X
and Y
from the tangent space at p
on the Manifold
M
using the Metric
G
. If G
is the default metric (see is_default_metric
) this is done using inner(M, p, X, Y)
, otherwise the local_metric
(M, p)
is employed as
where $G_p$ is the loal matrix representation of the Metric
G
.