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:

  1. to implement different metrics (e.g. in closed form) for one Manifold
  2. 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 β€” Type
Metric

Abstract type for the pseudo-Riemannian metric tensor $g$, a family of smoothly varying inner products on the tangent space. See inner.

source
Manifolds.MetricManifold β€” Type
MetricManifold{𝔽,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)

Generate the Manifold M as a manifold with the Metric G.

source
Manifolds.RiemannianMetric β€” Type
RiemannianMetric <: 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.

source

Implement Different Metrics on the same Manifold

In order to distinguish different metrics on one manifold, one can introduce two Metrics 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 β€” Method
exp(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.

source
Base.log β€” Method
log(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}.

source
Manifolds.christoffel_symbols_first β€” Method
christoffel_symbols_first(M::MetricManifold, p; backend=:default)

Compute the Christoffel symbols of the first kind in local coordinates. The Christoffel symbols are (in Einstein summation convention)

\[Ξ“_{ijk} = \frac{1}{2} \Bigl[g_{kj,i} + g_{ik,j} - g_{ij,k}\Bigr],\]

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)$.

source
Manifolds.christoffel_symbols_second β€” Method
christoffel_symbols_second(M::MetricManifold, x; backend=:default)

Compute the Christoffel symbols of the second kind in local coordinates. The Christoffel symbols are (in Einstein summation convention)

\[Ξ“^{l}_{ij} = g^{kl} Ξ“_{ijk},\]

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)$.

source
Manifolds.christoffel_symbols_second_jacobian β€” Method
christoffel_symbols_second_jacobian(M::MetricManifold, p; backend = :default)

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)$.

source
Manifolds.det_local_metric β€” Method
det_local_metric(M::MetricManifold, p)

Return the determinant of local matrix representation of the metric tensor $g$.

source
Manifolds.einstein_tensor β€” Method
einstein_tensor(M::MetricManifold, p; backend = :default)

Compute the Einstein tensor of the manifold M at the point p.

source
Manifolds.flat β€” Method
flat(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

\[X^β™­= G_p X,\]

where $G_p$ is the local matrix representation of G, see local_metric

source
Manifolds.gaussian_curvature β€” Method
gaussian_curvature(M::MetricManifold, x; backend = :default)

Compute the Gaussian curvature of the manifold M at the point x.

source
Manifolds.inverse_local_metric β€” Method
inverse_local_metric(M::MetricManifold, p)

Return the local matrix representation of the inverse metric (cometric) tensor, usually written $g^{ij}$.

source
Manifolds.is_default_metric β€” Method
is_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.

source
Manifolds.is_default_metric β€” Method
is_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.

source
Manifolds.local_metric β€” Method
local_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.

source
Manifolds.local_metric_jacobian β€” Method
local_metric_jacobian(M::MetricManifold, p; backend=:default)

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)$.

source
Manifolds.log_local_metric_density β€” Method
log_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}]|}$.

source
Manifolds.ricci_curvature β€” Method
ricci_curvature(M::MetricManifold, p; backend = :default)

Compute the Ricci scalar curvature of the manifold M at the point p.

source
Manifolds.ricci_tensor β€” Method
ricci_tensor(M::MetricManifold, p; backend = :default)

Compute the Ricci tensor, also known as the Ricci curvature tensor, of the manifold M at the point p.

source
Manifolds.riemann_tensor β€” Method
riemann_tensor(M::MetricManifold, p)

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)$.

source
Manifolds.sharp β€” Method
sharp(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

\[ΞΎ^β™― = G_p^{-1} ΞΎ,\]

where $G_p$ is the local matrix representation of G, i.e. one employs inverse_local_metric here to obtain $G_p^{-1}$.

source
ManifoldsBase.inner β€” Method
inner(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

\[g_p(X, Y) = ⟨X, G_p Y⟩,\]

where $G_p$ is the loal matrix representation of the Metric G.

source
Manifolds.solve_exp_ode β€” Method
solve_exp_ode(
    M::MetricManifold,
    p,
    X,
    tspan;
    backend = :default,
    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

\[\frac{d^2}{dt^2} p^k + Ξ“^k_{ij} \frac{d}{dt} p_i \frac{d}{dt} p_j = 0,\]

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.

Note

This function only works for Julia 1.1 or greater, when OrdinaryDiffEq.jl is loaded with

using OrdinaryDiffEq
source