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 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_indexcan 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_parametersconverts a point to its parameters with respect to the chart in a chart.get_pointconverts parameters (local coordinates) in a chart to the point that corresponds to them.induced_basisreturns a basis of a given vector space at a point induced by a chart $Ο$.transition_mapconverts 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 with charts that have values in the vector space π½βΏ for some value of n. π½ is a number system determined by an AbstractNumbers object.
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
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
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.
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.
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.
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.
Manifolds.christoffel_symbols_second β Method
christoffel_symbols_second(M::AbstractManifold, A::AbstractAtlas, i, a)Compute values of the Christoffel symbol of the second kind in chart i of atlas A at point with parameters a.
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
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
Manifolds.get_coordinates_induced_basis! β Method
get_coordinates_induced_basis!(M::AbstractManifold, c, p, X, B::InducedBasis{β, TangentSpaceType, <:AbstractAtlas};
backend::AbstractADType = AutoForwardDiff())Compute the coordinates of a tangent vector X at a point p on the manifold M in the induced basis B and store the result in c. This function uses automatic differentiation to compute the coordinates.
Arguments
M::AbstractManifold: The manifold on which the computation is performed.c: The output array where the coordinates of the tangent vector will be stored.p: The point on the manifold where the tangent vectorXis located.X: The tangent vector atpwhose coordinates are to be computed.B::InducedBasis{β, TangentSpaceType, <:AbstractAtlas}: The induced basis in which the coordinates are expressed.
Keyword arguments
backend::AbstractADType: The automatic differentiation backend used for computing derivatives (default:AutoForwardDiff()).
Returns
The result is stored in c, which contains the coordinates of the tangent vector X in the induced basis B.
Notes
- This function computes the coordinates by differentiating the chart map at the given point
pin the direction ofX. - The computation relies on automatic differentiation.
See also
Manifolds.get_default_atlas β Method
get_default_atlas(::AbstractManifold)Determine the default real-valued atlas for the given manifold.
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
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
Manifolds.get_vector_induced_basis! β Method
get_vector_induced_basis!(M::AbstractManifold, Y, p, Xc, B::InducedBasis{β, TangentSpaceType, <:AbstractAtlas};
backend::AbstractADType = AutoForwardDiff())Compute the tangent vector Y at a point p on the manifold M corresponding to the coordinates Xc in the induced basis B and store the result in Y. This function uses automatic differentiation to compute the tangent vector.
Arguments
M::AbstractManifold: The manifold on which the computation is performed.Y: The output tangent vector atpcorresponding to the coordinatesXcin the induced basis.p: The point on the manifold where the tangent vector is located.Xc: The coordinates of the tangent vector in the induced basisB.B::InducedBasis{β, TangentSpaceType, <:AbstractAtlas}: The induced basis in which the coordinatesXcare expressed.
Keyword arguments
backend::AbstractADType: The automatic differentiation backend used for computing derivatives (default:AutoForwardDiff()).
Returns
The result is stored in Y, which represents the tangent vector at p corresponding to the coordinates Xc in the induced basis B.
Notes
- This function computes the tangent vector by differentiating the chart map at the given point
pin the direction of the coordinatesXc. - The computation relies on automatic differentiation.
See also
Manifolds.induced_basis β Method
induced_basis(M::AbstractManifold, A::AbstractAtlas, i, 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 β Method
induced_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 β Method
inverse_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.inverse_local_metric β Method
inverse_local_metric(M::AbstractManifold, A::AbstractAtlas{β}, i, a)Compute the inverse of the local metric tensor on the manifold M at the point with parameters a in chart i of an AbstractAtlas A. The inverse local metric tensor is represented as a matrix, where each entry corresponds to the inverse of the inner product of basis vectors in the tangent space at the point with given parameters.
In contrast, inverse_local_metric(M::AbstractManifold, p, ::InducedBasis) requires passing a point instead of its parameters in a chart.
Arguments
M::AbstractManifold: The manifold on which the metric is computed.A::AbstractAtlas{β}: The atlas defining the charts and coordinate systems on the manifold.i: The index of the chart in the atlas.a: The parameters of the point in the chart.
Returns
A matrix representing the inverse of the local metric tensor at the point with given parameters.
See also
Manifolds.kretschmann_scalar β Method
kretschmann_scalar(M::AbstractManifold, A::AbstractAtlas, i, a; backend::AbstractADType = AutoForwardDiff())Compute the Kretschmann scalar $K = R_{abcd} R^{abcd}$ at the point given by coordinates a in chart i of atlas A on manifold M.
This implementation uses the Riemann tensor in the form $R^u_{ijk}$ (returned by riemann_tensor) and the inverse local metric g^{ij} (returned by inverse_local_metric) to form the full contraction:
\[ K = g^{u v} g^{i p} g^{j q} g^{k r} R^u_{i j k} R^v_{p q r}\]
Arguments
M::AbstractManifold: manifoldA::AbstractAtlas: atlas providing charts / induced basisi: chart index inAa: coordinates of the point in charti(lengthn)backend::AbstractADType: automatic-differentiation backend (defaultAutoForwardDiff())
Returns
Scalar (same element type as a) equal to the Kretschmann scalar at the point
Manifolds.levi_civita_affine_connection! β Method
levi_civita_affine_connection!(M::AbstractManifold, Zc, A::AbstractAtlas, i, a, Xc, Yc; backend::AbstractADType = AutoForwardDiff())Compute the Levi-Civita affine connection on the manifold M at a point with parameters a in chart i of an AbstractAtlas A. The connection is calculated for vectors with coefficients Xc and Yc in the induced basis, and the result is stored in Zc.
The Levi-Civita connection is computed using the metric tensor (inner called in a chart) of the manifold, ensuring that the connection is torsion-free and compatible with the metric. The computation involves the Christoffel symbols of the second kind, which are derived from the metric tensor and its derivatives using automatic differentiation. Note that this computation is relatively slow. Where performance matters, it should be replaced with custom-derived formulas.
Arguments
M::AbstractManifold: The manifold on which the Levi-Civita connection is computed.Zc: The output vector where the result of the connection is stored.A::AbstractAtlas: The atlas defining the charts and coordinate systems on the manifold.i: The index of the chart in the atlas.a: The parameters (local coordinates) of the point in the chart.Xc: The coefficients of the first vector in the induced basis.Yc: The coefficients of the second vector in the induced basis.
Keyword arguments
backend::AbstractADType: The automatic differentiation backend used for computing derivatives (default:AutoForwardDiff()).
Returns
The result is stored in Zc, which represents the Levi-Civita connection in the induced basis.
Notes
- The computation involves the inverse of the local metric tensor, which is used to raise indices.
- The directional derivatives of the metric tensor are computed using the specified automatic differentiation backend.
- The function assumes that the input vectors
XcandYcare expressed in the induced basis of the chart.
See also
Manifolds.local_metric β Method
local_metric(M::AbstractManifold, A::AbstractAtlas{β}, i, a)Compute the local metric tensor on the manifold M at the point with parameters a in chart i of an AbstractAtlas A. The local metric tensor is represented as a matrix, where each entry corresponds to the inner product of basis vectors in the tangent space at the point with given parameters.
In contrast, local_metric(M::AbstractManifold, p, ::InducedBasis) requires passing a point instead of its parameters in a chart.
Arguments
M::AbstractManifold: The manifold on which the metric is computed.A::AbstractAtlas{β}: The atlas defining the charts and coordinate systems on the manifold.i: The index of the chart in the atlas.a: The parameters of the point in the chart.
Returns
A matrix representing the local metric tensor at the point with given parameters.
See also
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.
Manifolds.ricci_curvature β Method
ricci_curvature(M::AbstractManifold, A::AbstractAtlas, i, a; backend=AutoForwardDiff())Compute the scalar Ricci curvature (Ricci scalar) of the manifold M at the point given by coordinates a in chart i of atlas A.
The scalar curvature is the trace of the Ricci tensor with respect to the inverse local metric:
\[ R = g^{ij} R_{ij}\]
math
Arguments
M::AbstractManifold: manifoldA::AbstractAtlas: atlas providing charts / induced basisi: chart index inAa: coordinates of the point in charti(lengthn)backend::AbstractADType: automatic-differentiation backend (defaultAutoForwardDiff())
Returns
- scalar (same element type as
a) equal to the Ricci scalar at the point
Manifolds.ricci_tensor β Method
ricci_tensor(M::AbstractManifold, A::AbstractAtlas, i, a; backend::AbstractADType=AutoForwardDiff())Compute the Ricci tensor of the manifold M at the point specified by coordinates a in chart i of atlas A.
The Ricci tensor is the contraction of the Riemann tensor:
\[ Ric_{p q} = R^u_{p u q}\]
Arguments
M::AbstractManifold: manifoldA::AbstractAtlas: atlas providing charts / induced basisi: chart index inAa: coordinates of the point in charti(lengthn)backend::AbstractADType: automatic-differentiation backend (defaultAutoForwardDiff())
Returns
nΓnmatrix with componentsRic[p, q](same element type asa)
See also
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
Manifolds.transition_map_diff! β Method
transition_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 β 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.
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.
ManifoldsBase.riemann_tensor β Method
riemann_tensor(M::AbstractManifold, A::AbstractAtlas, i, a;
backend::AbstractADType = AutoForwardDiff()Compute the Riemann curvature tensor of manifold M at the point given by parameters a in chart i of atlas A.
Returns a 4-dimensional array R of size (n,n,n,n) with components R[u,i,j,k] = R^u_{ijk}, where the first index is the contravariant (upper) index and the remaining three are covariant (lower) indices. The components satisfy, for coordinate vector fields e_i:
\[ R^u_{ijk} e_u = (β_{e_i} β_{e_j} - β_{e_j} β_{e_i} - β_{[e_i,e_j]}) e_k\]
Arguments
M::AbstractManifold: manifoldA::AbstractAtlas: atlas used for coordinates/induced basisi: chart index inAa: coordinates of the point in charti(lengthn)backend::AbstractADType: automatic-differentiation backend (defaultAutoForwardDiff())
Notes
- The default implementation computes connection coefficients via
affine_connectionand their directional derivatives (using automatic differentiation), so it can be expensive. Manifold-specific overrides yielding closed-form curvature are recommended for performance-critical code. - Use
riemann_tensor!(M, Wc, A, i, a, Xc, Yc, Zc)to compute the actionR(X,Y)Zon coordinate vectorsXc,Yc,Zcwithout constructing the full 4-tensor.
See also
ManifoldsBase.riemann_tensor β Method
riemann_tensor(M::AbstractManifold, A::AbstractAtlas, i, a, Xc, Yc, Zc;
backend::AbstractADType = AutoForwardDiff())Compute the action of the Riemann curvature tensor R on tangent vectors with coordinates Xc, Yc and Zc at the point specified by parameters a in chart i of atlas A on manifold M.
This function returns the vector W (in induced-chart coordinates) given by(R(X, Y) Z), i.e. the result of applying the curvature operator toZc`.
Arguments
M::AbstractManifold: manifoldA::AbstractAtlas: atlas providing charts / induced basisi: chart index inAa: coordinates of the point in charti(lengthn)Xc, Yc, Zc: coordinates of tangent vectors X, Y, Z in the chart-induced basis
Keyword arguments
backend::AbstractADType = AutoForwardDiff(): AD backend used when numerical derivatives are required
Returns
- A vector (same shape/type as
Xc) containing coordinates of $(R(X, Y) Z)$ in the induced basis.
Notes
- The default implementation builds the full 4-index Riemann tensor using the chart affine connection and its derivatives, then contracts with
Xc,Yc,Zc. This is convenient but can be expensive; prefer manifold-specific overrides that compute R(X,Y)Z directly for performance-critical code. - Inputs are expected to be expressed in the chart-induced basis associated with
Aandi. - The function uses automatic differentiation (via
backend) when computing directional derivatives of connection coefficients.
See also
riemann_tensor(M, A, i, a)which returns the full 4-tensor- []
affine_connection](@ref) used to obtain connection coefficients (seelevi_civita_affine_connection!for a generic implementation)
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 β Type
RieszRepresenterCotangentVector(M::AbstractManifold, p, X)Cotangent vector in Riesz representer form on manifold M at point p with Riesz representer X.
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$
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$
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.
Manifolds.estimate_distance_from_bvp β Function
estimate_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 β Function
solve_chart_exp_ode(
M::AbstractManifold, a, Xc, A::AbstractAtlas, i0;
solver=AutoVern9(Rodas5()),
final_time::Real=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. The geodesic is solved up to time final_time (by default equal to 1).
Chart switching
If the solution exceeds the domain of chart i0 (which is detected using the check_chart_switch function with additional keyword arguments check_chart_switch_kwargs), a new chart is selected using get_chart_index on the final point in the old chart.
Returned value
The function returns an object of type StitchedChartSolution{:Exp} to represent the geodesic.
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.
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.