Introduction

Riemannian Manifolds

All manifolds inherit from Manifold to store their main properties, which is most prominently the manifold dimension and the name of the manifold. This will be extended in the future, for example properties denoting whether the manifold is explicitly given in the sense of a closed form exponential and logarithmic map for example, or only approximately.

A Riemannian manifold in Manopt.jl consist of three types:

Manopt.ManifoldType.

An abstract manifold $\mathcal M$ to keep global information on a specific manifold

source
Manopt.MPointType.

An abstract point $x$ on a Manifold $\mathcal M$.

source
Manopt.TVectorType.

A tangent vector $\xi \in T_x\mathcal M$ at a MPoint point $x$ on a Manifold $\mathcal M$.

source

List of available Manifolds

Furthermore there are two types accompanying each manifold – a point on the manifold inheriting from MPoint and the tangent vector TVector. For both the term manifold is shortened to M for concise naming. Each manifold also inherits such a short abbreviation, see Abbr. in the following table.

Manifold $\mathcal M$FileAbbr.Comment
A manifold $\mathcal M$Manifold.jlM
$1$-sphere $\mathbb S^1$Circle.jlS1represented as angles $x\in[-\pi,\pi)$
Euclidean space $\mathbb R^n$Euclidean.jlRn$n$-dimensional Euclidean space $\mathbb R^n$
Grassmannian manifold $\mathrm{Gr}(k,n)$Grassmannian.jlGrembedded in $\mathbb R^{n\times k}$
$n$-dim. Hyperbolic space $\mathbb H^n$Hyperbolic.jlHnembedded in $\mathbb R^{n+1}$
special orthogonal group $\mathrm{SO}(n)$Rotations.jlSOrepresented as rotation matrices
$n$-sphere $\mathbb S^n$Sphere.jlSnembedded in $\mathbb R^{n+1}$
Stiefel $\mathrm{St}(k,n)$Stiefel.jlStcontains both the real- ad the complex-valued case
symmetric matrices $\mathcal{Sym}(n)$Symmetric.jlSym$n\times n$ symmetric matrices
symmetric positive definite matrices $\mathcal P(n)$SymmetricPositiveDefinite.jlSPD$n\times n$ symmetric positive matrices using the affine metric

If you're missing your favorite manifold, give us a note on Github.

Special Types of Manifolds

Manifolds build upon Manifolds

Manifold $\mathcal M$FileAbbr.Comment
Power manifoldPower.jlPowBuilds $\mathcal N^{d_1\times\cdot\times d_k}$ of any manifold $\mathcal N$
Product manifoldProduct.jlProdBuild the product manifold $\mathcal N_1\times\cdots\times\mathcal N_k$ of manifolds
Tangent bundleTangentBundle.jlTBtangent bundle of a manifold, i.e. the set of all tuples $(x,\xi), \xi \in T_x\mathcal M$, $x\in\mathcal M$ with the induced metric.

for more details see Combined Manifolds

Special Properties of Manifolds

Special types of manifolds are introduced by SimpleTraits.jl. They can be used to clarify that a manifold possesses a certain property. For example two points on a matrix manifold can be multiplied, though the result is not necessarily a point on the manifold anymore. Traits have to goals here: Provide functions that are common for all manifolds of such a type (e.g. the for Lie groups) as a common interface and to specify certain functions or solvers for these certain types, that for example take advantage of then.

Embedded Manifold

IsEmbeddedM{X}

An abstract Manifold that is embedded in some Euclidean space. These manifolds may have projections and converters for gradient and Hessian.

source
IsEmbeddedP{X}

An abstract MPoint belonging to an embedded manifold.

source
IsEmbeddedV{X}

An abstract TVector belonging to an embedded manifold.

source

Lie Group Manifold

IsLieGroupM{X}

Indicates that X is a Manifold with a Lie group structure. This also introdcues a group operation of two MPoints of X.

source
IsLieGroupP{X}

An abstract MPoint belonging to a Lie group manifold.

source
IsLieGroupV{X}

An abstract TVector belonging to a Lie group manifold.

source
Manopt.:⊗Function.
⊗(x,y)

the binary operator x ⊗ y represents the Lie group action for a Lie group Manifold, see IsLieGroupM.

source

Matrix Manifold

IsMatrixM{X}

An abstract Manifold to represent a manifold whose points are matrices. For these manifolds the usual operators (+,-,*) are overloaded for points. Furthermore, the transpose is also overloaded, though it returns the matrix, since the dimensions mit be different for rectangular matrices.

source
IsMatrixP{X}

An abstract Manifold Point belonging to a matrix manifold.

source
IsMatrixTV{X}

An abstract Manifold Point belonging to a matrix manifold.

source

Functions that need to be implemented for a Manifold

If you plan to implement a new manifold within Manopt.jl, the following functions should be implemented. If you only implement a few of these functions, not all algorithms might work. all these functions have a fallback providing an error message if the function is not (yet) implemented. Otherwise, for example, if the field of the inner representant of MPoint or TVector is the field .value of your data structure, the default implementation of getValue directly works. In the following list M <: Manifold the manifold type represents the manifold Q,P <: MPoint the type of a point on the new manifold, T <: TVector a corresponding tangent vector in a suitable tangent space,

Manopt.addNoiseMethod.
addNoise(M,x)

add noise to a MPoint x on the Manifold M by using the randomTVector method and doing an exponential step. Optional parameters, like the type of noise and parameters for the noise may be given and are just passed on-

source
Manopt.distanceMethod.
distance(M,x,y)

computes the gedoesic distance between two MPoints x and y on a Manifold M.

source
LinearAlgebra.dotMethod.
dot(M, x, ξ, ν)

Computes the inner product of two TVectors ξ and ν from the tangent space at the MPoint x on the Manifold M.

source
Base.expMethod.
exp(M,x,ξ,[t=1.0])

computes the exponential map at an MPoint x for the TVector ξ on the Manifold M. The optional parameter t can be used to shorten ξ to .

source
Manopt.getValueMethod.
getValue(x)

get the value representing the MPoint x. This function defaults to returning x.value; if your representation is different, you should implement this function for your type

source
Manopt.getValueMethod.
getValue(ξ)

get the value representing the TVector ξ. This function defaults to returning ξ.value; if your representation is different, you should implement this function for your type

source
Base.logMethod.
log(M,x,y)

computes the TVector in the tangent space $T_x\mathcal M$ at the MPoint x such that the corresponding geodesic reaches the MPoint y after time 1 on the Manifold M.

source
manifoldDimension(x)

return the dimension of the manifold M the point x belongs to.

source
manifoldDimension(M)

returns the dimension of the manifold M.

source
LinearAlgebra.normMethod.
norm(M,x,ξ)

computes the length of a TVector ξ in the tangent space of the MPoint x on the Manifold M induced by the inner product.

source
parallelTransport(M,x,y,ξ)

Parallel transport of a vector ξ given at the tangent space $T_x\mathcal M$ of x to the tangent space $T_y\mathcal M$ at y along the geodesic form x to y. If the geodesic is not unique, this function takes the same choice as geodesic.

source
randomMPoint(M,x [,:Gaussian,options...])

generate a random MPoint on the Manifold M by falling back to the default :Gaussian noise with the default standard deviation on the specific manifold.

source
randomTVector(M,x [,:Gaussian,options...])

generate a random tangent vector at MPoint x on the Manifold M using :Gaussian noise where options usually contain the standard deviation σ on the specific manifold.

source
Manopt.tangentONBMethod.
(Ξ,κ) = tangentONB(M,x,y)

compute an ONB within the tangent space $T_x\mathcal M$ such that $\xi=\log_xy$ is the first vector and compute the eigenvalues of the curvature tensor $R(\Xi,\dot g)\dot g$, where $g=g_{x,\xi}$ is the geodesic with $g(0)=x$, $\dot g(0) = \xi$, i.e. $\kappa_1$ corresponding to $\Xi_1=\xi$ is zero.

See also

jacobiField, adjointJacobiField.

source
Manopt.tangentONBMethod.
(Ξ,κ) = tangentONB(M,x,ξ)

compute an ONB within the tangent space $T_x\mathcal M$ such that $\xi$ is the first vector and compute the eigenvalues of the curvature tensor $R(\Xi,\dot g)\dot g$, where $g=g_{x,\xi}$ is the geodesic with $g(0)=x$, $\dot g(0) = \xi$, i.e. $\kappa_1$ corresponding to $\Xi_1=\xi$ is zero.

See also

jacobiField, adjointJacobiField.

source
typeofMPoint(ξ)

return the MPoint belonging to the TVector type of ξ.

source
typeofMPoint(T)

return the MPoint belonging to the TVector type T.

source
typeofTVector(x)

return the TVector belonging to the MPoint type of x.

source
typeofTVector(P)

returns the TVector belonging to the MPoint type P.

source
typicalDistance(M)

returns the typical distance on the Manifold M, which is for example the longest distance in a unit cell or injectivity radius.

source
Manopt.zeroTVectorMethod.
ξ = zeroTVector(M,x)

returns a zero vector in the tangent space $T_x\mathcal M$ of the MPoint $x\in\mathcal M$ on the Manifold M.

source

Functions implemented for a general manifold

the following base functions are implemented for general manifolds and are based on the functions from the last section

ζ = adjointJacobiField(M,x,y,t,η,w)

compute the AdjointJacobiField $J$ along the geodesic $g_{x,y}$ on the manifold $\mathcal M$ with initial conditions (depending on the application) $\eta\in T_{g(t;x,y)\mathcal M}$ and weights $\beta$. The result is a vector $\zeta \in T_x\mathcal M$. The main difference to jacobiField is the, that the input $\eta$ and the output $\zeta$ switched tangent spaces.

For detais see jacobiField

source
Manopt.geodesicFunction.
geodesic(M,x,y)

return a function to evaluate the geodesic connecting the two MPoints x and y on the Manifold M.

source
geodesic(M,x,y,n)

return vector containing the equispaced n sample-values along the geodesic connecting the two MPoints x and y on the Manifold M.

source
geodesic(M,x,y,t)

return the point along the geodesic from MPoint x to y given by at value t (in [0,1]) on the Manifold M

source
geodesic(M,x,y,T)

return vector containing the MPoint along the geodesic from MPoint x to y on the Manifold M specified by the points from the vector T (of numbers between 0 and 1).

source
Manopt.jacobiFieldFunction.
ζ = jacobiField(M,x,y,t,η,β)

compute the jacobiField $J$ along the geodesic $g_{x,y}$ on the Manifold M $\mathcal M$ with initial conditions (depending on the application) $\eta\in T_x\mathcal M$ and weights $\beta$. The result is a TVector in $\zeta \in T_{g(t;x,y)}\mathcal M$.

See also

adjointJacobiField

source
Statistics.meanFunction.
y = mean(M,x)

compute the Riemannian center of mass of the data given by the vector of MPoints x on the Manifold M. calculates the Riemannian Center of Mass (Karcher mean) of the input data x as an Array of MPoints on the Manifold with a steepestDescent or a cyclicProximalPoint.

Optional

  • initialValue – (x[1]) the value to initialize the algorithm to
  • method – (:GradientDescent) symbol indicating the algorithm to use, so the second variant is :CyclicProximalPoint
  • weights – (1/n) compute a weighted Karcher mean, i.e. the dault is to set all weights to be 1/n where n is the length of x.

as well as optional parameters that are passed down to the corresponding algorithm

source
Statistics.medianFunction.
y = median(M,x)

compute the median of the data given by the vector of MPoints x on the Manifold M. calculates the Riemannian Center of Mass (Karcher mean) of the input data x as an Array of MPoints on the Manifold with a cyclicProximalPoint.

Optional

  • initialValue – (x[1]) the value to initialize the algorithm to
  • method – (:CyclicProximalPoint) symbol indicating the algorithm to use
  • weights – (1/n) compute a weighted Karcher mean, i.e. the dault is to set all weights to be 1/n where n is the length of x.

as well as optional parameters that are passed down to the corresponding algorithm

source
Manopt.midPointFunction.
midPoint(M,x,y,z)

compute the mid point between x and y. If there is more than one mid point of (not neccessarily miniizing) geodesics (i.e. on the sphere), the one nearest to z.

source
midPoint(M,x,y)

compute the (geodesic) mid point of the two MPoints x and y on the Manifold M. If the geodesic is not unique, either a deterministic choice is returned or an error is raised. For the deteministic choixe, see midPoint(M,x,y,z), the mid point closest to a third MPoint z.

source
Manopt.reflectionFunction.
y = reflection(M,p,x)

reflect the MPoint x at MPoint p, i.e. compute $y = R_p(x) = \exp_p(-\log_px)$. On Euclidean space this results in the point reflection $R_p(x) = p - (x-p) = 2p-x$.

Arguments

  • M – a Manifold $\mathcal M$
  • p – an MPoint $p\in\mathcal M$ to relfect at
  • x – an MPoint $x\in\mathcal M$ that is reflected

Output

  • y – the resulting reflection.
source
Manopt.varianceFunction.
variance(M,x)

returns the variance of the vector x of MPoints on the Manifold M

source

A decorator for checks and validation

In order to check and/or validate tangent vectors, the decorator pattern TVectorE is available for any subtype of TVector as follows

Manopt.TVectorEType.
TVectorE <: MPoint

A decorator pattern based extension of TVector to additionally store the base point. The decorator is then used to verify, that exp and dot are only called with correct base points.

Constructors

TVectorE(ξ,x [,v=true])

constructs an extended tangential vector based on the TVector ξ with base MPoint x with optional validation v. If none of the first two arguments is an extended element, v defaults to true, otherwise, the one which is an extended is inherited or the && of both validations.

source

together with a small helper MPointE that indicates to the log to returns an extended tangent vector as soon as one of its arguments is an extended manifold point.

Manopt.MPointEType.
MPointE <: MPoint

A decorator pattern based extension of MPoint to identify when to switch to the extended TVectorE for functions just working on points, e.g. log. The constructor avoids multiple encapsualtions of extensions.

Constructors

MPointE(x [,v=true])

the point can constructed by extending an existing MPoint. optionally, the validation can be turned off (but is trueby default). If x is a MPointE the default of v is taken from x.

source

for these two data items, the following additional features that are activated

Inheritance

Basic functions like exp, log, parallelTransport, and randomTVector return an extended MPointE or TVectorE whenever one of its arguments is an extended input. This enables, that setting (only) one input for a calculation to an extended version, this property propagates this into all the algorithm.

Note that this might increase memory usage and hence reduce performance, since for any TVectorE internally stores both a TVector as well as its base MPoint (as extension).

Checks

Manopt.checkBasePointFunction.
checkBasePoint(ξ,ν)

checks, whether the base of two tangent vectors is identical, if both tangent vectors are of type TVectorE. If one of them is not an extended vector, the function returns true, expecting the tangent vector implicitly to be correct.

source
checkBasePoint(ξ,x)

checks, whether the base of the tangent vector ξ is x. If ξ is not an extended tangent vector TVectorE the function returns true, assuming the base implicitly to be correct

source

For extended data decorators, whenever possible in the basic functions listed above a checkBasePoint completely automatically performed. For example, when calling exp(M,x,ξ), as soon as ξ is an extended vector, checkBasePoint(x,ξ) is called before performing the original exponential map.

This way, as many checks are performed, whether corresponding points and vectors or two vectors involved have the correct base points.

Validation

Every extended type carries a further boolean validation, whose default is true, i.e. to perform validation. activating validation one needs to implement the following two functions, otherwise, a lot of warnings might occur.

validateMPoint(M,x)

check, whether the data in the MPoint x is a valid point on the Manifold M. This is used to validate parameters and results during computations using MPointEs. Note that the default fallback is just a warning that no validation is available.

The function should throw an error if x is not point on the manifold M, otherwise it should return true.

source
validateTVector(M,x,ξ)

check, whether the data in the TVector ξ is a valid tangent vector #to the MPoint x on the Manifold M. This is used to validate parameters and results during computations when using MPointEs. Note that the default fallback is just a warning that no validation is available.

Available validations should throw an error if x is not on M or ξ is not in the tangent space of x. If ξ is valid, the function returns true.

source

whenever possible (see checkBasePoint above). Since such a validation might not be available for your favorite manifold, you can deactivate validation by setting the boolean to false. Every new extended type inherits the false, whenever one of its part (i.e. either the TVector or MPoint) has validation set to false.

So while checking as often as possible, this feature can easily be deactivated.

Internals

To access the inner value of MPointE and the base stored in TVectorE you can use

Base.stripFunction.
strip(ξ)

returns the internal TVector of an extended tangent vector TVectorE.

source
strip(x)

returns the MPoint the MPointE x stores internally. If applied to an already non-extended MPoint, nothing happens.

source
Manopt.getBasePointFunction.
getBasePoint(ξ)

returns the base point of an extended tangent vector. To continue promotion of the extended type, the result is always a MPointE. To eliminate the decorator, use strip.

source

Furthermore the following functions are mapping to the internally stored data and encapsulate the results with the extended variant if applicable

as well as mathematical operators on tangent vectors and comparison operators.