Atlases and charts

Atlases on an $n$-dimensional manifold $mathcal M$are collections of charts $\mathcal A = \{(U_i, Ο†_i) \colon i \in I\}$, where $I$ is a (finite or infinte) index family, such that $U_i \subseteq \mathcal M$ is an open set and each chart $Ο†_i: U_i β†’ ℝ^n$ is a homeomorphism. This means, that $Ο†_i$ is bijective – sometimes also called one-to-one and onto - and continuous, and its inverse $Ο†_i^{-1}$ is continuous as well. The inverse $Ο†_i^{-1}$ is called (local) parametrization. The resulting parameters $a=Ο†(p)$ of $p$ (with respect to the chart $Ο†$) are in the literature also called β€œ(local) coordinates”. To distinguish the parameter $a$ from get_coordinates in a basis, we use the terminology parameter in this package.

For an atlas $\mathcal A$ we further require that

\[\displaystyle\bigcup_{i\in I} U_i = \mathcal M.\]

We say that $Ο†_i$ is a chart about $p$, if $p\in U_i$. An atlas provides a connection between a manifold and the Euclidean space $ℝ^n$, since locally, a chart about $p$ can be used to identify its neighborhood (as long as you stay in $U_i$) with a subset of a Euclidean space. Most manifolds we consider are smooth, i.e. any change of charts $Ο†_i \circ Ο†_j^{-1}: ℝ^n β†’ ℝ^n$, where $i,j\in I$, is a smooth function. These changes of charts are also called transition maps.

Most operations on manifolds in Manifolds.jl avoid operating in a chart through appropriate embeddings and formulas derived for particular manifolds, though atlases provide the most general way of working with manifolds. Compared to these approaches, using an atlas is often more technical and time-consuming. They are extensively used in metric-related functions on MetricManifolds.

Atlases are represented by objects of subtypes of AbstractAtlas. There are no type restrictions for indices of charts in atlases.

Operations using atlases and charts are available through the following functions:

  • get_chart_index can be used to select an appropriate chart for the neighborhood of a given point $p$. This function should work deterministically, i.e. for a fixed $p$ always return the same chart.
  • get_parameters converts a point to its parameters with respect to the chart in a chart.
  • get_point converts parameters (local coordinates) in a chart to the point that corresponds to them.
  • induced_basis returns a basis of a given vector space at a point induced by a chart $Ο†$.
  • transition_map converts coordinates of a point between two charts, e.g. computes $Ο†_i\circ Ο†_j^{-1}: ℝ^n β†’ ℝ^n$, $i,j\in I$.

While an atlas could store charts as explicit functions, it is favourable, that the [get_parameters] actually implements a chart $Ο†$, get_point its inverse, the prametrization $Ο†^{-1}$.

Manifolds.AbstractAtlas β€” Type
AbstractAtlas{𝔽}

An abstract class for atlases whith charts that have values in the vector space 𝔽ⁿ for some value of n. 𝔽 is a number system determined by an AbstractNumbers object.

source
Manifolds.InducedBasis β€” Type
InducedBasis(vs::VectorSpaceType, A::AbstractAtlas, i)

The basis induced by chart with index i from an AbstractAtlas A of vector space of type vs.

For the vs a TangentSpace this works as follows:

Let $n$ denote the dimension of the manifold $\mathcal M$.

Let the parameter $a=Ο†_i(p) ∈ \mathbb R^n$ and $j∈\{1,…,n\}$. We can look at the $j$th parameter curve $b_j(t) = a + te_j$, where $e_j$ denotes the $j$th unit vector. Using the parametrisation we obtain a curve $c_j(t) = Ο†_i^{-1}(b_j(t))$ which fulfills $c(0) = p$.

Now taking the derivative(s) with respect to $t$ (and evaluate at $t=0$), we obtain a tangent vector for each $j$ corresponding to an equivalence class of curves (having the same derivative) as

\[X_j = [c_j] = \frac{\mathrm{d}}{\mathrm{d}t} c_i(t) \Bigl|_{t=0}\]

and the set $\{X_1,\ldots,X_n\}$ is the chart-induced basis of $T_p\mathcal M$.

See also

VectorSpaceType, AbstractBasis

source
Manifolds.RetractionAtlas β€” Type
RetractionAtlas{
    𝔽,
    TRetr<:AbstractRetractionMethod,
    TInvRetr<:AbstractInverseRetractionMethod,
    TBasis<:AbstractBasis,
} <: AbstractAtlas{𝔽}

An atlas indexed by points on a manifold, $\mathcal M = I$ and parameters (local coordinates) are given in $T_p\mathcal M$. This means that a chart $Ο†_p = \mathrm{cord}\circ\mathrm{retr}_p^{-1}$ is only locally defined (around $p$), where $\mathrm{cord}$ is the decomposition of the tangent vector into coordinates with respect to the given basis of the tangent space, cf. get_coordinates. The parametrization is given by $Ο†_p^{-1}=\mathrm{retr}_p\circ\mathrm{vec}$, where $\mathrm{vec}$ turns the basis coordinates into a tangent vector, cf. get_vector.

In short: The coordinates with respect to a basis are used together with a retraction as a parametrization.

See also

AbstractAtlas, AbstractInverseRetractionMethod, AbstractRetractionMethod, AbstractBasis

source
LinearAlgebra.norm β€” Method
norm(M::AbstractManifold, A::AbstractAtlas, i, a, Xc)

Calculate norm on manifold M at point with parameters a in chart i of an AbstractAtlas A of vector with coefficients Xc in induced basis.

source
Manifolds.affine_connection! β€” Method
affine_connection!(M::AbstractManifold, Zc, A::AbstractAtlas, i, a, Xc, Yc)

Calculate affine connection on manifold M at point with parameters a in chart i of an an AbstractAtlas A of vectors with coefficients Zc and Yc in induced basis and save the result in Zc.

source
Manifolds.affine_connection β€” Method
affine_connection(M::AbstractManifold, A::AbstractAtlas, i, a, Xc, Yc)

Calculate affine connection on manifold M at point with parameters a in chart i of AbstractAtlas A of vectors with coefficients Xc and Yc in induced basis.

source
Manifolds.check_chart_switch β€” Method
check_chart_switch(M::AbstractManifold, A::AbstractAtlas, i, a)

Determine whether chart should be switched when an operation in chart i from an AbstractAtlas A reaches parameters a in that chart.

By default false is returned.

source
Manifolds.get_chart_index β€” Method
get_chart_index(M::AbstractManifold, A::AbstractAtlas, i, a)

Select a chart from an AbstractAtlas A for manifold M that is suitable for representing the neighborhood of point with parametrization a in chart i. This selection should be deterministic, although different charts may be selected for arbitrarily close but distinct points.

See also

get_default_atlas

source
Manifolds.get_chart_index β€” Method
get_chart_index(M::AbstractManifold, A::AbstractAtlas, p)

Select a chart from an AbstractAtlas A for manifold M that is suitable for representing the neighborhood of point p. This selection should be deterministic, although different charts may be selected for arbitrarily close but distinct points.

See also

get_default_atlas

source
Manifolds.get_parameters β€” Method
get_parameters(M::AbstractManifold, A::AbstractAtlas, i, p)

Calculate parameters (local coordinates) of point p on manifold M in chart from an AbstractAtlas A at index i. This function is hence an implementation of the chart $Ο†_i(p), i\in I$. The parameters are in the number system determined by A. If the point $p\notin U_i$ is not in the domain of the chart, this method should throw an error.

See also

get_point, get_chart_index

source
Manifolds.get_point β€” Method
get_point(M::AbstractManifold, A::AbstractAtlas, i, a)

Calculate point at parameters (local coordinates) a on manifold M in chart from an AbstractAtlas A at index i. This function is hence an implementation of the inverse $Ο†_i^{-1}(a), i\in I$ of a chart, also called a parametrization.

See also

get_parameters, get_chart_index

source
Manifolds.local_metric β€” Method
local_metric(M::AbstractManifold, p, B::InducedBasis)

Compute the local metric tensor for vectors expressed in terms of coordinates in basis B on manifold M. The point p is not checked.

source
Manifolds.transition_map β€” Method
transition_map(M::AbstractManifold, A_from::AbstractAtlas, i_from, A_to::AbstractAtlas, i_to, a)
transition_map(M::AbstractManifold, A::AbstractAtlas, i_from, i_to, a)

Given coordinates a in chart (A_from, i_from) of a point on manifold M, returns coordinates of that point in chart (A_to, i_to). If A_from and A_to are equal, A_to can be omitted.

Mathematically this function is the transition map or change of charts, but it might even be between two atlases $A_{\text{from}} = \{(U_i,Ο†_i)\}_{i\in I}$ and $A_{\text{to}} = \{(V_j,\psi_j)\}_{j\in J}$, and hence $I, J$ are their index sets. We have $i_{\text{from}}\in I$, $i_{\text{to}}\in J$.

This method then computes

\[\bigl(\psi_{i_{\text{to}}}\circ Ο†_{i_{\text{from}}}^{-1}\bigr)(a)\]

Note that, similarly to get_parameters, this method should fail the same way if $V_{i_{\text{to}}}\cap U_{i_{\text{from}}}=\emptyset$.

See also

AbstractAtlas, get_parameters, get_point

source
Manifolds.transition_map_diff β€” Method
transition_map_diff(M::AbstractManifold, A::AbstractAtlas, i_from, a, c, i_to)

Compute differential of transition map from chart i_from to chart i_to from an AbstractAtlas A on manifold M at point with parameters a on tangent vector with coordinates c in the induced basis.

source
ManifoldsBase.inner β€” Method
inner(M::AbstractManifold, A::AbstractAtlas, i, a, Xc, Yc)

Calculate inner product on manifold M at point with parameters a in chart i of an atlas A of vectors with coefficients Xc and Yc in induced basis.

source

Cotangent space and musical isomorphisms

Related to atlases, there is also support for the cotangent space and coefficients of cotangent vectors in bases of the cotangent space.

Functions sharp and flat implement musical isomorphisms for arbitrary vector bundles.

Manifolds.flat β€” Method
flat(M::AbstractManifold, p, X)

Compute the flat isomorphism (one of the musical isomorphisms) of tangent vector X from the vector space of type M at point p from the underlying AbstractManifold.

The function can be used for example to transform vectors from the tangent bundle to vectors from the cotangent bundle $β™­ : T\mathcal M β†’ T^{*}\mathcal M$

source
Manifolds.sharp β€” Method
sharp(M::AbstractManifold, p, ΞΎ)

Compute the sharp isomorphism (one of the musical isomorphisms) of vector ΞΎ from the vector space M at point p from the underlying AbstractManifold.

The function can be used for example to transform vectors from the cotangent bundle to vectors from the tangent bundle $β™― : T^{*}\mathcal M β†’ T\mathcal M$

source

Computations in charts

Manifolds.IntegratorTerminatorNearChartBoundary β€” Type
IntegratorTerminatorNearChartBoundary{TKwargs}

An object for determining the point at which integration of a differential equation in a chart on a manifold should be terminated for the purpose of switching a chart.

The value stored in check_chart_switch_kwargs will be passed as keyword arguments to check_chart_switch. By default an empty tuple is stored.

source
Manifolds.solve_chart_exp_ode β€” Function
solve_chart_exp_ode(
    M::AbstractManifold,
    a,
    Xc,
    A::AbstractAtlas,
    i0;
    solver=AutoVern9(Rodas5()),
    final_time=1.0,
    check_chart_switch_kwargs=NamedTuple(),
    kwargs...,
)

Solve geodesic ODE on a manifold M from point of coordinates a in chart i0 from an AbstractAtlas A in direction of coordinates Xc in the induced basis.

source
Manifolds.solve_chart_log_bvp β€” Function
solve_chart_log_bvp(
    M::AbstractManifold,
    a1,
    a2,
    A::AbstractAtlas,
    i;
    solver=MIRK4(),
    dt::Real=0.05,
    kwargs...,
)

Solve the BVP corresponding to geodesic calculation on AbstractManifold M, between points with parameters a1 and a2 in a chart i of an AbstractAtlas A using solver solver. Geodesic Ξ³ is sampled at time interval dt, with Ξ³(0) = a1 and Ξ³(1) = a2.

source
Manifolds.solve_chart_parallel_transport_ode β€” Function
solve_chart_parallel_transport_ode(
    M::AbstractManifold,
    a,
    Xc,
    A::AbstractAtlas,
    i0,
    Yc;
    solver=AutoVern9(Rodas5()),
    check_chart_switch_kwargs=NamedTuple(),
    final_time=1.0,
    kwargs...,
)

Parallel transport vector with coordinates Yc along geodesic on a manifold M from point of coordinates a in a chart i0 from an AbstractAtlas A in direction of coordinates Xc in the induced basis.

source