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 infinite) 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 MetricManifold
s.
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
β TypeAbstractAtlas{π½}
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.
Manifolds.InducedBasis
β TypeInducedBasis(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
Manifolds.RetractionAtlas
β TypeRetractionAtlas{
π½,
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
LinearAlgebra.norm
β Methodnorm(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.
Manifolds.affine_connection!
β Methodaffine_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
.
Manifolds.affine_connection
β Methodaffine_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.
Manifolds.check_chart_switch
β Methodcheck_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.
Manifolds.get_chart_index
β Methodget_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
Manifolds.get_chart_index
β Methodget_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
Manifolds.get_default_atlas
β Methodget_default_atlas(::AbstractManifold)
Determine the default real-valued atlas for the given manifold.
Manifolds.get_parameters
β Methodget_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
Manifolds.get_point
β Methodget_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
Manifolds.induced_basis
β Methodinduced_basis(M::AbstractManifold, A::AbstractAtlas, i, p, VST::VectorSpaceType)
Basis of vector space of type VST
at point p
from manifold M
induced by chart (A
, i
).
See also
Manifolds.induced_basis
β Methodinduced_basis(::AbstractManifold, A::AbstractAtlas, i, VST::VectorSpaceType = TangentSpaceType())
Get the basis induced by chart with index i
from an AbstractAtlas
A
of vector space of type vs
. Returns an object of type InducedBasis
.
See also
Manifolds.inverse_chart_injectivity_radius
β Methodinverse_chart_injectivity_radius(M::AbstractManifold, A::AbstractAtlas, i)
Injectivity radius of get_point
for chart i
from an AbstractAtlas
A
of a manifold M
.
Manifolds.local_metric
β Methodlocal_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.
Manifolds.transition_map
β Methodtransition_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
Manifolds.transition_map_diff!
β Methodtransition_map_diff!(M::AbstractManifold, c_out, A::AbstractAtlas, i_from, a, c, i_to)
Compute transition_map_diff
on given arguments and save the result in c_out
.
Manifolds.transition_map_diff
β Methodtransition_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.
ManifoldsBase.inner
β Methodinner(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.
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.RieszRepresenterCotangentVector
β TypeRieszRepresenterCotangentVector(M::AbstractManifold, p, X)
Cotangent vector in Riesz representer form on manifold M
at point p
with Riesz representer X
.
Manifolds.flat
β Methodflat(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$
Manifolds.sharp
β Methodsharp(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$
Computations in charts
Manifolds.IntegratorTerminatorNearChartBoundary
β TypeIntegratorTerminatorNearChartBoundary{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.
Manifolds.estimate_distance_from_bvp
β Functionestimate_distance_from_bvp(
M::AbstractManifold,
a1,
a2,
A::AbstractAtlas,
i;
solver=MIRK4(),
dt=0.05,
kwargs...,
)
Estimate distance between points on AbstractManifold
M with parameters a1
and a2
in chart i
of AbstractAtlas
A
using solver solver
, employing solve_chart_log_bvp
to solve the geodesic BVP.
Manifolds.solve_chart_exp_ode
β Functionsolve_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.
Manifolds.solve_chart_log_bvp
β Functionsolve_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.
Manifolds.solve_chart_parallel_transport_ode
β Functionsolve_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.