Functions on manifolds

This page collects several basic functions on manifolds.

The exponential map, the logarithmic map, and geodesics

Geodesics are the generalizations of a straight line to manifolds, i.e. their intrinsic acceleration is zero. Together with geodesics one also obtains the exponential map and its inverse, the logarithmic map. Informally speaking, the exponential map takes a vector (think of a direction and a length) at one point and returns another point, which lies towards this direction at distance of the specified length. The logarithmic map does the inverse, i.e. given two points, it tells which vector “points towards” the other point.

Base.expMethod
exp(M::AbstractManifold, p, X)
exp(M::AbstractManifold, p, X, t::Number = 1)

Compute the exponential map of tangent vector X, optionally scaled by t, at point p from the manifold AbstractManifold M, i.e.

\[\exp_p X = γ_{p,X}(1),\]

where $γ_{p,X}$ is the unique geodesic starting in $γ(0)=p$ such that $\dot γ(0) = X$.

See also shortest_geodesic, retract.

source
Base.logMethod
log(M::AbstractManifold, p, q)

Compute the logarithmic map of point q at base point p on the AbstractManifold M. The logarithmic map is the inverse of the exponential map. Note that the logarithmic map might not be globally defined.

See also inverse_retract.

source
ManifoldsBase.exp!Method
exp!(M::AbstractManifold, q, p, X)
exp!(M::AbstractManifold, q, p, X, t::Number = 1)

Compute the exponential map of tangent vector X, optionally scaled by t, at point p from the manifold AbstractManifold M. The result is saved to q.

If you want to implement exponential map for your manifold, you should implement the method with t, that is exp!(M::MyManifold, q, p, X, t::Number).

See also exp.

source
ManifoldsBase.geodesic!Method
geodesic!(M::AbstractManifold, Q, p, X, T::AbstractVector) -> AbstractVector

Get the geodesic with initial point p and velocity X on the AbstractManifold M. A geodesic is a curve of zero acceleration. That is for the curve $γ_{p,X}: I → \mathcal M$, with $γ_{p,X}(0) = p$ and $\dot γ_{p,X}(0) = X$ a geodesic further fulfills

\[∇_{\dot γ_{p,X}(t)} \dot γ_{p,X}(t) = 0,\]

i.e. the curve is acceleration free with respect to the Riemannian metric. This function evaluates the geodeic at time points t fom T in place of Q.

source
ManifoldsBase.geodesic!Method
geodesic!(M::AbstractManifold, q, p, X, t::Real)

Get the geodesic with initial point p and velocity X on the AbstractManifold M. A geodesic is a curve of zero acceleration. That is for the curve $γ_{p,X}: I → \mathcal M$, with $γ_{p,X}(0) = p$ and $\dot γ_{p,X}(0) = X$ a geodesic further fulfills

\[∇_{\dot γ_{p,X}(t)} \dot γ_{p,X}(t) = 0,\]

i.e. the curve is acceleration free with respect to the Riemannian metric. This function evaluates the geodeic at t in place of q.

source
ManifoldsBase.geodesic!Method
geodesic!(M::AbstractManifold, p, X) -> Function

Get the geodesic with initial point p and velocity X on the AbstractManifold M. A geodesic is a curve of zero acceleration. That is for the curve $γ_{p,X}: I → \mathcal M$, with $γ_{p,X}(0) = p$ and $\dot γ_{p,X}(0) = X$ a geodesic further fulfills

\[∇_{\dot γ_{p,X}(t)} \dot γ_{p,X}(t) = 0,\]

i.e. the curve is acceleration free with respect to the Riemannian metric. This yields that the curve has constant velocity and is locally distance-minimizing.

This function returns a function (q,t) of (time) t that mutates q`.

source
ManifoldsBase.geodesicMethod
geodesic(M::AbstractManifold, p, X, T::AbstractVector) -> AbstractVector

Evaluate the geodesic $γ_{p,X}: I → \mathcal M$, with $γ_{p,X}(0) = p$ and $\dot γ_{p,X}(0) = X$ a geodesic further fulfills

\[∇_{\dot γ_{p,X}(t)} \dot γ_{p,X}(t) = 0,\]

at time points t from T.

source
ManifoldsBase.geodesicMethod
geodesic(M::AbstractManifold, p, X, t::Real)

Evaluate the geodesic $γ_{p,X}: I → \mathcal M$, with $γ_{p,X}(0) = p$ and $\dot γ_{p,X}(0) = X$ a geodesic further fulfills

\[∇_{\dot γ_{p,X}(t)} \dot γ_{p,X}(t) = 0,\]

at time t.

source
ManifoldsBase.geodesicMethod
geodesic(M::AbstractManifold, p, X) -> Function

Get the geodesic with initial point p and velocity X on the AbstractManifold M. A geodesic is a curve of zero acceleration. That is for the curve $γ_{p,X}: I → \mathcal M$, with $γ_{p,X}(0) = p$ and $\dot γ_{p,X}(0) = X$ a geodesic further fulfills

\[∇_{\dot γ_{p,X}(t)} \dot γ_{p,X}(t) = 0,\]

i.e. the curve is acceleration free with respect to the Riemannian metric. This yields, that the curve has constant velocity that is locally distance-minimizing.

This function returns a function of (time) t.

source
ManifoldsBase.log!Method
log!(M::AbstractManifold, X, p, q)

Compute the logarithmic map of point q at base point p on the AbstractManifold M. The result is saved to X. The logarithmic map is the inverse of the exp!onential map. Note that the logarithmic map might not be globally defined.

see also log and inverse_retract!,

source
ManifoldsBase.shortest_geodesic!Method
shortest_geodesic!(M::AbstractManifold, R, p, q, T::AbstractVector) -> AbstractVector

Evaluate a geodesic $γ_{p,q}(t)$ whose length is the shortest path between the points pand q, where $γ_{p,q}(0)=p$ and $γ_{p,q}(1)=q$ at all t from T in place of R. When there are multiple shortest geodesics, a deterministic choice will be taken.

source
ManifoldsBase.shortest_geodesic!Method
shortest_geodesic!(M::AabstractManifold, r, p, q, t::Real)

Evaluate a geodesic $γ_{p,q}(t)$ whose length is the shortest path between the points pand q, where $γ_{p,q}(0)=p$ and $γ_{p,q}(1)=q$ at t in place of r. When there are multiple shortest geodesics, a deterministic choice will be taken.

source
ManifoldsBase.shortest_geodesic!Method
shortest_geodesic!(M::AbstractManifold, p, q) -> Function

Get a geodesic $γ_{p,q}(t)$ whose length is the shortest path between the points pand q, where $γ_{p,q}(0)=p$ and $γ_{p,q}(1)=q$. When there are multiple shortest geodesics, a deterministic choice will be returned.

This function returns a function (r,t) -> ... of time t which works in place of r.

Further variants

shortest_geodesic!(M::AabstractManifold, r, p, q, t::Real)
shortest_geodesic!(M::AbstractManifold, R, p, q, T::AbstractVector) -> AbstractVector

mutate (and return) the point r and the vector of points R, respectively, returning the point at time t or points at times t in T along the shortest geodesic.

source
ManifoldsBase.shortest_geodesicMethod
shortest_geodesic(M::AbstractManifold, p, q, T::AbstractVector) -> AbstractVector

Evaluate a geodesic $γ_{p,q}(t)$ whose length is the shortest path between the points pand q, where $γ_{p,q}(0)=p$ and $γ_{p,q}(1)=q$ at time points T. When there are multiple shortest geodesics, a deterministic choice will be returned.

source
ManifoldsBase.shortest_geodesicMethod
shortest_geodesic(M::AabstractManifold, p, q, t::Real)

Evaluate a geodesic $γ_{p,q}(t)$ whose length is the shortest path between the points pand q, where $γ_{p,q}(0)=p$ and $γ_{p,q}(1)=q$ at time t. When there are multiple shortest geodesics, a deterministic choice will be returned.

source
ManifoldsBase.shortest_geodesicMethod
shortest_geodesic(M::AbstractManifold, p, q) -> Function

Get a geodesic $γ_{p,q}(t)$ whose length is the shortest path between the points pand q, where $γ_{p,q}(0)=p$ and $γ_{p,q}(1)=q$. When there are multiple shortest geodesics, a deterministic choice will be returned.

This function returns a function of time, which may be a Real or an AbstractVector.

source

Parallel transport

While moving vectors from one base point to another is the identity in the Euclidean space – or in other words all tangent spaces (directions one can “walk” into) are the same. This is different on a manifold.

If we have two points $p,q ∈ \mathcal M$, we take a $c: [0,1] → \mathcal M$ connecting the two points, i.e. $c(0) = p$ and $c(1) = q$. this could be a (or the) geodesic. If we further consider a vector field $X: [0,1] → T\mathcal M$, i.e. where $X(t) ∈ T_{c(t)}\mathcal M$. Then the vector field is called parallel if its covariant derivative $\frac{\mathrm{D}}{\mathrm{d}t}X(t) = 0$ for all $t∈ |0,1]$.

If we now impose a value for $X=X(0) ∈ T_p\mathcal M$, we obtain an ODE with an initial condition. The resulting value $X(1) ∈ T_q\mathcal M$ is called the parallel transport of X along $c$ or in case of a geodesic the _parallel transport of X from p to q.

ManifoldsBase.parallel_transport_alongMethod
Y = parallel_transport_along(M::AbstractManifold, p, X, c)

Compute the parallel transport of the vector X from the tangent space at p along the curve c.

To be precise let $c(t)$ be a curve $c(0)=p$ for vector_transport_along $\mathcal P^cY$

THen In the result $Y\in T_p\mathcal M$ is the vector $X$ from the tangent space at $p=c(0)$ to the tangent space at $c(1)$.

Let $Z\colon [0,1] \to T\mathcal M$, $Z(t)\in T_{c(t)}\mathcal M$ be a smooth vector field along the curve $c$ with $Z(0) = Y$, such that $Z$ is parallel, i.e. its covariant derivative $\frac{\mathrm{D}}{\mathrm{d}t}Z$ is zero. Note that such a $Z$ always exists and is unique.

Then the parallel transport is given by $Z(1)$.

source

Further functions on manifolds

General functions provided by the interface

Base.angleMethod
angle(M::AbstractManifold, p, X, Y)

Compute the angle between tangent vectors X and Y at point p from the AbstractManifold M with respect to the inner product from inner.

source
Base.copyto!Method
copyto!(M::AbstractManifold, Y, p, X)

Copy the value(s) from X to Y, where both are tangent vectors from the tangent space at p on the AbstractManifold M. This function defaults to calling copyto!(Y, X), but it might be useful to overwrite the function at the level, where also information from p and M can be accessed.

source
Base.copyto!Method
copyto!(M::AbstractManifold, q, p)

Copy the value(s) from p to q, where both are points on the AbstractManifold M. This function defaults to calling copyto!(q, p), but it might be useful to overwrite the function at the level, where also information from M can be accessed.

source
Base.isapproxMethod
isapprox(M::AbstractManifold, p, X, Y; error:Symbol=:none; kwargs...)

Check if vectors X and Y tangent at p from AbstractManifold M are approximately equal.

The optional positional argument can be used to get more information for the case that the result is false, if the concrete manifold provides such information. Currently the following are supported

  • :error - throws an error if isapprox evaluates to false, providing possibly a more detailed error. Note that this turns isapprox basically to an @assert.
  • :info – prints the information in an @info
  • :warn – prints the information in an @warn
  • :none (default) – the function just returns true/false

By default these informations are collected by calling check_approx.

Keyword arguments can be used to specify tolerances.

source
Base.isapproxMethod
isapprox(M::AbstractManifold, p, q; error::Symbol=:none, kwargs...)

Check if points p and q from AbstractManifold M are approximately equal.

The keyword argument can be used to get more information for the case that the result is false, if the concrete manifold provides such information. Currently the following are supported

  • :error - throws an error if isapprox evaluates to false, providing possibly a more detailed error. Note that this turns isapprox basically to an @assert.
  • :info – prints the information in an @info
  • :warn – prints the information in an @warn
  • :none (default) – the function just returns true/false

Keyword arguments can be used to specify tolerances.

source
Base.randMethod
Random.rand(M::AbstractManifold, [d::Integer]; vector_at=nothing)
Random.rand(rng::AbstractRNG, M::AbstractManifold, [d::Integer]; vector_at=nothing)

Generate a random point on manifold M (when vector_at is nothing) or a tangent vector at point vector_at (when it is not nothing).

Optionally a random number generator rng to be used can be specified. An optional integer d indicates that a vector of d points or tangent vectors is to be generated.

Note

Usually a uniform distribution should be expected for compact manifolds and a Gaussian-like distribution for non-compact manifolds and tangent vectors, although it is not guaranteed. The distribution may change between releases.

rand methods for specific manifolds may take additional keyword arguments.

source
ManifoldsBase.WeingartenMethod
Weingarten(M, p, X, V)

Compute the Weingarten map $\mathcal W_p\colon T_p\mathcal M × N_p\mathcal M \to T_p\mathcal M$, where $N_p\mathcal M$ is the orthogonal complement of the tangent space $T_p\mathcal M$ of the embedded submanifold $\mathcal M$, where we denote the embedding by $\mathcal E$.

The Weingarten map can be defined by restricting the differential of the orthogonal projection $\operatorname{proj}_{T_p\mathcal M}\colon T_p \mathcal E \to T_p\mathcal M$ with respect to the base point $p$, i.e. defining

\[\mathcal P_X := D_p\operatorname{proj}_{T_p\mathcal M}(Y)[X], \qquad Y \in T_p \mathcal E, X \in T_p\mathcal M,\]

the Weingarten map can be written as $\mathcal W_p(X,V) = \mathcal P_X(V)$.

The Weingarten map is named after Julius Weingarten (1836–1910).

source
ManifoldsBase.allocateMethod
allocate(a)
allocate(a, dims::Integer...)
allocate(a, dims::Tuple)
allocate(a, T::Type)
allocate(a, T::Type, dims::Integer...)
allocate(a, T::Type, dims::Tuple)
allocate(M::AbstractManifold, a)
allocate(M::AbstractManifold, a, dims::Integer...)
allocate(M::AbstractManifold, a, dims::Tuple)
allocate(M::AbstractManifold, a, T::Type)
allocate(M::AbstractManifold, a, T::Type, dims::Integer...)
allocate(M::AbstractManifold, a, T::Type, dims::Tuple)

Allocate an object similar to a. It is similar to function similar, although instead of working only on the outermost layer of a nested structure, it maps recursively through outer layers and calls similar on the innermost array-like object only. Type T is the new number element type number_eltype, if it is not given the element type of a is retained. The dims argument can be given for non-nested allocation and is forwarded to the function similar.

It's behavior can be overridden by a specific manifold, for example power manifold with nested replacing representation can decide that allocate for Array{<:SArray} returns another Array{<:SArray} instead of Array{<:MArray}, as would be done by default.

source
ManifoldsBase.allocate_onMethod
allocate_on(M::AbstractManifold, [T:::Type])
allocate_on(M::AbstractManifold, F::FiberType, [T:::Type])

Allocate a new point on manifold M with optional type given by T. Note that T is not number element type as in allocate but rather the type of the entire point to be returned.

If F is provided, then an element of the corresponding fiber is allocated, assuming it is independent of the base point.

To allocate a tangent vector, use ``

Example

julia> using ManifoldsBase

julia> M = ManifoldsBase.DefaultManifold(4)
DefaultManifold(4; field = ℝ)

julia> allocate_on(M)
4-element Vector{Float64}:
 0.0
 0.0
 0.0
 0.0

julia> allocate_on(M, Array{Float64})
4-element Vector{Float64}:
 0.0
 0.0
 0.0
 0.0

julia> allocate_on(M, TangentSpaceType())
4-element Vector{Float64}:
 0.0
 0.0
 0.0
 0.0

julia> allocate_on(M, TangentSpaceType(), Array{Float64})
4-element Vector{Float64}:
 0.0
 0.0
 0.0
 0.0
source
ManifoldsBase.base_manifoldFunction
base_manifold(M::AbstractManifold, depth = Val(-1))

Return the internally stored AbstractManifold for decorated manifold M and the base manifold for vector bundles or power manifolds. The optional parameter depth can be used to remove only the first depth many decorators and return the AbstractManifold from that level, whether its decorated or not. Any negative value deactivates this depth limit.

source
ManifoldsBase.default_typeMethod
default_type(M::AbstractManifold, ft::FiberType)

Get the default type of points from the fiber ft of the fiber bundle based on manifold M. For example, call default_type(MyManifold(), TangentSpaceType()) to get the default type of a tangent vector.

source
ManifoldsBase.distanceMethod
distance(M::AbstractManifold, p, q)

Shortest distance between the points p and q on the AbstractManifold M, i.e.

\[d(p,q) = \inf_{γ} L(γ),\]

where the infimum is over all piecewise smooth curves $γ: [a,b] \to \mathcal M$ connecting $γ(a)=p$ and $γ(b)=q$ and

\[L(γ) = \displaystyle\int_{a}^{b} \lVert \dotγ(t)\rVert_{γ(t)} \mathrm{d}t\]

is the length of the curve $γ$.

If $\mathcal M$ is not connected, i.e. consists of several disjoint components, the distance between two points from different components should be $∞$.

source
ManifoldsBase.embed!Method
embed!(M::AbstractManifold, Y, p, X)

Embed a tangent vector X at a point p on the AbstractManifold M into the ambient space and return the result in Y. This method is only available for manifolds where implicitly an embedding or ambient space is given. Additionally, embed! includes changing data representation, if applicable, i.e. if the tangents on M are not represented in the same way as tangents on the embedding, the representation is changed accordingly. This is the case for example for Lie groups, when tangent vectors are represented in the Lie algebra. The embedded tangents are then in the tangent spaces of the embedded base points.

The default is set in such a way that it assumes that the points on M are represented in their embedding (for example like the unit vectors in a space to represent the sphere) and hence embedding also for tangent vectors is the identity by default.

See also: EmbeddedManifold, project!

source
ManifoldsBase.embed!Method
embed!(M::AbstractManifold, q, p)

Embed point p from the AbstractManifold M into an ambient space. This method is only available for manifolds where implicitly an embedding or ambient space is given. Not implementing this function means, there is no proper embedding for your manifold. Additionally, embed might include changing data representation, if applicable, i.e. if points on M are not represented in the same way as their counterparts in the embedding, the representation is changed accordingly.

The default is set in such a way that it assumes that the points on M are represented in their embedding (for example like the unit vectors in a space to represent the sphere) and hence embedding in the identity by default.

If you have more than one embedding, see EmbeddedManifold for defining a second embedding. If your point p is already represented in some embedding, see AbstractDecoratorManifold how you can avoid reimplementing code from the embedded manifold

See also: EmbeddedManifold, project!

source
ManifoldsBase.embedMethod
embed(M::AbstractManifold, p, X)

Embed a tangent vector X at a point p on the AbstractManifold M into an ambient space. This method is only available for manifolds where implicitly an embedding or ambient space is given. Not implementing this function means, there is no proper embedding for your tangent space(s).

Additionally, embed might include changing data representation, if applicable, i.e. if tangent vectors on M are not represented in the same way as their counterparts in the embedding, the representation is changed accordingly.

The default is set in such a way that memory is allocated and embed!(M, Y, p. X) is called.

If you have more than one embedding, see EmbeddedManifold for defining a second embedding. If your tangent vector X is already represented in some embedding, see AbstractDecoratorManifold how you can avoid reimplementing code from the embedded manifold

See also: EmbeddedManifold, project

source
ManifoldsBase.embedMethod
embed(M::AbstractManifold, p)

Embed point p from the AbstractManifold M into the ambient space. This method is only available for manifolds where implicitly an embedding or ambient space is given. Additionally, embed includes changing data representation, if applicable, i.e. if the points on M are not represented in the same way as points on the embedding, the representation is changed accordingly.

The default is set in such a way that memory is allocated and embed!(M, q, p) is called.

See also: EmbeddedManifold, project

source
ManifoldsBase.embed_projectMethod
embed_project(M::AbstractManifold, p, X)

Embed vector X tangent at p from manifold M an project it back to tangent space at p. For points from that tangent space this is identity but in case embedding is defined for tangent vectors from outside of it, this can serve as a way to for example remove numerical inaccuracies caused by some algorithms.

source
ManifoldsBase.embed_projectMethod
embed_project(M::AbstractManifold, p)

Embed p from manifold M an project it back to M. For points from M this is identity but in case embedding is defined for points outside of M, this can serve as a way to for example remove numerical inaccuracies caused by some algorithms.

source
ManifoldsBase.injectivity_radiusMethod
injectivity_radius(M::AbstractManifold)

Infimum of the injectivity radii injectivity_radius(M,p) of all points p on the AbstractManifold.

injectivity_radius(M::AbstractManifold, p)

Return the distance $d$ such that exp(M, p, X) is injective for all tangent vectors shorter than $d$ (i.e. has an inverse).

injectivity_radius(M::AbstractManifold[, x], method::AbstractRetractionMethod)
injectivity_radius(M::AbstractManifold, x, method::AbstractRetractionMethod)

Distance $d$ such that retract(M, p, X, method) is injective for all tangent vectors shorter than $d$ (i.e. has an inverse) for point p if provided or all manifold points otherwise.

In order to dispatch on different retraction methods, please either implement _injectivity_radius(M[, p], m::T) for your retraction R or specifically injectivity_radius_exp(M[, p]) for the exponential map. By default the variant with a point p assumes that the default (without p) can ve called as a lower bound.

source
ManifoldsBase.is_pointMethod
is_point(M::AbstractManifold, p; error::Symbol = :none, kwargs...)
is_point(M::AbstractManifold, p, throw_error::Bool; kwargs...)

Return whether p is a valid point on the AbstractManifold M. By default the function calls check_point, which returns an ErrorException or nothing.

How to report a potential error can be set using the error= keyword

  • :error - throws an error if p is not a point
  • :info - displays the error message as an @info
  • :warn - displays the error message as a @warning
  • :none (default) – the function just returns true/false

all other symbols are equivalent to error=:none.

The second signature is a shorthand, where the boolean is used for error=:error (true) and error=:none (default, false). This case ignores the error= keyword

source
ManifoldsBase.is_vectorMethod
is_vector(M::AbstractManifold, p, X, check_base_point::Bool=true; error::Symbol=:none, kwargs...)
is_vector(M::AbstractManifold, p, X, check_base_point::Bool=true, throw_error::Boolean; kwargs...)

Return whether X is a valid tangent vector at point p on the AbstractManifold M. Returns either true or false.

If check_base_point is set to true, this function also (first) calls is_point on p. Then, the function calls check_vector and checks whether the returned value is nothing or an error.

How to report a potential error can be set using the error= keyword

  • :error - throws an error if X is not a tangent vector and/or p is not point

^ :info - displays the error message as an @info

  • :warn - displays the error message as a @warning.
  • :none - (default) the function just returns true/false

all other symbols are equivalent to error=:none

The second signature is a shorthand, where throw_error is used for error=:error (true) and error=:none (default, false). This case ignores the error= keyword.

source
ManifoldsBase.mid_point!Method
mid_point!(M::AbstractManifold, q, p1, p2)

Calculate the middle between the two point p1 and p2 from manifold M. By default uses log, divides the vector by 2 and uses exp!. Saves the result in q.

source
ManifoldsBase.mid_pointMethod
mid_point(M::AbstractManifold, p1, p2)

Calculate the middle between the two point p1 and p2 from manifold M. By default uses log, divides the vector by 2 and uses exp.

source
ManifoldsBase.riemann_tensorMethod
riemann_tensor(M::AbstractManifold, p, X, Y, Z)

Compute the value of the Riemann tensor $R(X_f,Y_f)Z_f$ at point p, where $X_f$, $Y_f$ and $Z_f$ are vector fields defined by parallel transport of, respectively, X, Y and Z to the desired point. All computations are performed using the connection associated to manifold M.

The formula reads $R(X_f,Y_f)Z_f = \nabla_X\nabla_Y Z - \nabla_Y\nabla_X Z - \nabla_{[X, Y]}Z$, where $[X, Y]$ is the Lie bracket of vector fields.

Note that some authors define this quantity with inverse sign.

source
ManifoldsBase.sectional_curvatureMethod
sectional_curvature(M::AbstractManifold, p, X, Y)

Compute the sectional curvature of a manifold $\mathcal M$ at a point $p \in \mathcal M$ on two linearly independent tangent vectors at $p$. The formula reads

\[ \kappa_p(X, Y) = \frac{⟨R(X, Y, Y), X⟩_p}{\lVert X \rVert^2_p \lVert Y \rVert^2_p - ⟨X, Y⟩^2_p} \]

where $R(X, Y, Y)$ is the riemann_tensor on $\mathcal M$.

Input

  • M: a manifold $\mathcal M$
  • p: a point $p \in \mathcal M$
  • X: a tangent vector $X \in T_p \mathcal M$
  • Y: a tangent vector $Y \in T_p \mathcal M$
source
ManifoldsBase.sectional_curvature_maxMethod
sectional_curvature_max(M::AbstractManifold)

Upper bound on sectional curvature of manifold M. The formula reads

\[\omega = \operatorname{sup}_{p\in\mathcal M, X\in T_p\mathcal M, Y\in T_p\mathcal M, ⟨X, Y⟩ ≠ 0} \kappa_p(X, Y)\]

source
ManifoldsBase.sectional_curvature_minMethod
sectional_curvature_min(M::AbstractManifold)

Lower bound on sectional curvature of manifold M. The formula reads

\[\omega = \operatorname{inf}_{p\in\mathcal M, X\in T_p\mathcal M, Y\in T_p\mathcal M, ⟨X, Y⟩ ≠ 0} \kappa_p(X, Y)\]

source
ManifoldsBase.zero_vector!Method
zero_vector!(M::AbstractManifold, X, p)

Save to X the tangent vector from the tangent space $T_p\mathcal M$ at p that represents the zero vector, i.e. such that retracting X to the AbstractManifold M at p produces p.

source
ManifoldsBase.zero_vectorMethod
zero_vector(M::AbstractManifold, p)

Return the tangent vector from the tangent space $T_p\mathcal M$ at p on the AbstractManifold M, that represents the zero vector, i.e. such that a retraction at p produces p.

source

Internal functions

While you should always add your documentation to functions from the last section, some of the functions dispatch onto functions on layer III. These are the ones you usually implement for your manifold – unless there is no lower level function called, like for the manifold_dimension.

Base.convertMethod
convert(T::Type, M::AbstractManifold, p, X)

Convert vector X tangent at point p from manifold M to type T.

source
Base.convertMethod
convert(T::Type, M::AbstractManifold, p)

Convert point p from manifold M to type T.

source
ManifoldsBase._isapproxMethod
_isapprox(M::AbstractManifold, p, X, Y; kwargs...)

An internal function for testing whether tangent vectors X and Y from tangent space at point p from manifold M are approximately equal. Returns either true or false and does not support errors like isapprox.

For more details see documentation of check_approx.

source
ManifoldsBase._isapproxMethod
_isapprox(M::AbstractManifold, p, q; kwargs...)

An internal function for testing whether points p and q from manifold M are approximately equal. Returns either true or false and does not support errors like isapprox.

For more details see documentation of check_approx.

source
ManifoldsBase._pick_basic_allocation_argumentMethod
_pick_basic_allocation_argument(::AbstractManifold, f, x...)

Pick which one of elements of x should be used as a basis for allocation in the allocate_result(M::AbstractManifold, f, x...) method. This can be specialized to, for example, skip Identity arguments in Manifolds.jl group-related functions.

source
ManifoldsBase.allocate_resultMethod
allocate_result(M::AbstractManifold, f, x...)

Allocate an array for the result of function f on AbstractManifold M and arguments x... for implementing the non-modifying operation using the modifying operation.

Usefulness of passing a function is demonstrated by methods that allocate results of musical isomorphisms.

source
ManifoldsBase.check_approxMethod
check_approx(M::AbstractManifold, p, q; kwargs...)
check_approx(M::AbstractManifold, p, X, Y; kwargs...)

Check whether two elements are approximately equal, either p, q on the AbstractManifold or the two tangent vectors X, Y in the tangent space at p are approximately the same. The keyword arguments kwargs can be used to set tolerances, similar to Julia's isapprox.

This function might use isapprox from Julia internally and is similar to isapprox, with the difference that is returns an ApproximatelyError if the two elements are not approximately equal, containting a more detailed description/reason. If the two elements are approximalely equal, this method returns nothing.

This method is an internal function and is called by isapprox whenever the user specifies an error= keyword therein. _isapprox is another related internal function. It is supposed to provide a fast true/false decision whether points or vectors are equal or not, while check_approx also provides a textual explanation. If no additional explanation is needed, a manifold may just implement a method of _isapprox, while it should also implement check_approx if a more detailed explanation could be helpful.

source
ManifoldsBase.check_pointMethod
check_point(M::AbstractManifold, p; kwargs...) -> Union{Nothing,String}

Return nothing when p is a point on the AbstractManifold M. Otherwise, return an error with description why the point does not belong to manifold M.

By default, check_point returns nothing, i.e. if no checks are implemented, the assumption is to be optimistic for a point not deriving from the AbstractManifoldPoint type.

source
ManifoldsBase.check_sizeMethod
check_size(M::AbstractManifold, p)
check_size(M::AbstractManifold, p, X)

Check whether p has the right representation_size for a AbstractManifold M. Additionally if a tangent vector is given, both p and X are checked to be of corresponding correct representation sizes for points and tangent vectors on M.

By default, check_size returns nothing, i.e. if no checks are implemented, the assumption is to be optimistic.

source
ManifoldsBase.check_vectorMethod
check_vector(M::AbstractManifold, p, X; kwargs...) -> Union{Nothing,String}

Check whether X is a valid tangent vector in the tangent space of p on the AbstractManifold M. An implementation does not have to validate the point p. If it is not a tangent vector, an error string should be returned.

By default, check_vector returns nothing, i.e. if no checks are implemented, the assumption is to be optimistic for tangent vectors not deriving from the TVector type.

source

Approximation Methods

ManifoldsBase.ExtrinsicEstimationType
ExtrinsicEstimation{T} <: AbstractApproximationMethod

Method for estimation in the ambient space with a method of type T and projecting the result back to the manifold.

source
ManifoldsBase.GeodesicInterpolationWithinRadiusType
GeodesicInterpolationWithinRadius{T} <: AbstractApproximationMethod

Method for estimation based on geodesic interpolation that is restricted to some radius

Constructor

GeodesicInterpolationWithinRadius(radius::Real)
source
ManifoldsBase.default_approximation_methodMethod
default_approximation_method(M::AbstractManifold, f)
default_approximation_method(M::AbtractManifold, f, T)

Specify a default estimation method for an AbstractManifold and a specific function f and optionally as well a type T to distinguish different (point or vector) representations on M.

By default, all functions f call the signature for just a manifold. The exceptional functions are:

source

Error Messages

This interface introduces a small set of own error messages.

ManifoldsBase.ApproximatelyErrorType
ApproximatelyError{V,S} <: Exception

Store an error that occurs when two data structures, e.g. points or tangent vectors.

Fields

  • val amount the two approximate elements are apart – is set to NaN if this is not known
  • msg a message providing more detail about the performed test and why it failed.

Constructors

ApproximatelyError(val::V, msg::S) where {V,S}

Generate an Error with value val and message msg.

ApproximatelyError(msg::S) where {S}

Generate a message without a value (using val=NaN internally) and message msg.

source
ManifoldsBase.ComponentManifoldErrorType
CompnentError{I,E} <: Exception

Store an error that occured in a component, where the additional index is stored.

Fields

  • index::I index where the error occured`
  • error::E error that occured.
source
ManifoldsBase.ManifoldDomainErrorType
ManifoldDomainError{<:Exception} <: Exception

An error to represent a nested (Domain) error on a manifold, for example if a point or tangent vector is invalid because its representation in some embedding is already invalid.

source