Differentiation
Documentation for Manifolds.jl's methods and types for finite differences and automatic differentiation.
Differentiation backends
Manifolds.AbstractDiffBackend — TypeAbstractDiffBackendAn abstract type for diff backends. See FiniteDifferencesBackend for an example.
Manifolds.CurrentDiffBackend — TypeCurrentDiffBackend(backend::AbstractDiffBackend)A mutable struct for storing the current differentiation backend in a global constant _current_default_differential_backend.
See also
AbstractDiffBackend, default_differential_backend, set_default_differential_backend!
Manifolds._derivative — Function_derivative(f, t[, backend::AbstractDiffBackend])Compute the derivative of a callable f at time t computed using the given backend, an object of type Manifolds.AbstractDiffBackend. If the backend is not explicitly specified, it is obtained using the function default_differential_backend.
This function calculates plain Euclidean derivatives, for Riemannian differentiation see for example differential.
Not specifying the backend explicitly will usually result in a type instability and decreased performance.
Manifolds._gradient — Function_gradient(f, p[, backend::AbstractDiffBackend])Compute the gradient of a callable f at point p computed using the given backend, an object of type AbstractDiffBackend. If the backend is not explicitly specified, it is obtained using the function default_differential_backend.
This function calculates plain Euclidean gradients, for Riemannian gradient calculation see for example gradient.
Not specifying the backend explicitly will usually result in a type instability and decreased performance.
Manifolds._hessian — Function_hessian(f, p[, backend::AbstractDiffBackend])Compute the Hessian of a callable f at point p computed using the given backend, an object of type AbstractDiffBackend. If the backend is not explicitly specified, it is obtained using the function default_differential_backend.
This function calculates plain Euclidean Hessian.
Not specifying the backend explicitly will usually result in a type instability and decreased performance.
Manifolds._jacobian — Function_jacobian(f, p[, backend::AbstractDiffBackend])Compute the Jacobian of a callable f at point p computed using the given backend, an object of type AbstractDiffBackend. If the backend is not explicitly specified, it is obtained using the function default_differential_backend.
This function calculates plain Euclidean Jacobians, for Riemannian Jacobian calculation see for example gradient.
Not specifying the backend explicitly will usually result in a type instability and decreased performance.
Manifolds.default_differential_backend — Methoddefault_differential_backend() -> AbstractDiffBackendGet the default differentiation backend.
Manifolds.set_default_differential_backend! — Methodset_default_differential_backend!(backend::AbstractDiffBackend)Set current backend for differentiation to backend.
Manifolds._current_default_differential_backend — Constant_current_default_differential_backendThe instance of Manifolds.CurrentDiffBackend that stores the globally default differentiation backend.
Further differentiation backends and features are available in ManifoldDiff.jl.
FiniteDifferenes.jl
Manifolds.FiniteDifferencesBackend — TypeFiniteDifferencesBackend(method::FiniteDifferenceMethod = central_fdm(5, 1))Differentiation backend based on the FiniteDifferences package.
Riemannian differentiation backends
Manifolds.AbstractRiemannianDiffBackend — TypeAbstractRiemannianDiffBackendAn abstract type for backends for differentiation.
Manifolds.RiemannianProjectionBackend — TypeRiemannianProjectionBackend <: AbstractRiemannianDiffBackendThis backend computes the differentiation in the embedding, which is currently limited to the gradient. Let $mathcal M$ denote a manifold embedded in some $R^m$, where $m$ is usually (much) larger than the manifold dimension. Then we require three tools
- A function $f̃: ℝ^m → ℝ$ such that its restriction to the manifold yields the cost function $f$ of interest.
- A
projectfunction to project tangent vectors from the embedding (at $T_pℝ^m$) back onto the tangent space $T_p\mathcal M$. This also includes possible changes of the representation of the tangent vector (e.g. in the Lie algebra or in a different data format). - A
change_representerfor non-isometrically embedded manifolds, i.e. where the tangent space $T_p\mathcal M$ of the manifold does not inherit the inner product from restriction of the inner product from the tangent space $T_pℝ^m$ of the embedding
For more details see [AbsilMahonySepulchre2008], Section 3.6.1 for a derivation on submanifolds.
Manifolds.TangentDiffBackend — TypeTangentDiffBackend <: AbstractRiemannianDiffBackendA backend that uses tangent spaces and bases therein to derive an intrinsic differentiation scheme.
Since it works in tangent spaces at argument and function value, methods might require a retraction and an inverse retraction as well as a basis.
In the tangent space itself, this backend then employs an (Euclidean) AbstractDiffBackend
Constructor
TangentDiffBackend(diff_backend)where diff_backend is an AbstractDiffBackend to be used on the tangent space.
With the keyword arguments
retractionan AbstractRetractionMethod (ExponentialRetractionby default)inverse_retractionan AbstractInverseRetractionMethodLogarithmicInverseRetractionby default)basis_argan AbstractBasis (DefaultOrthogonalBasisby default)basis_valan AbstractBasis (DefaultOrthogonalBasisby default)
Manifolds.differential — Methoddifferential(M::AbstractManifold, f, t::Real, backend::AbstractDiffBackend)
differential!(M::AbstractManifold, f, X, t::Real, backend::AbstractDiffBackend)Compute the Riemannian differential of a curve $f: ℝ\to M$ on a manifold M represented by function f at time t using the given backend. It is calculated as the tangent vector equal to $\mathrm{d}f_t(t)[1]$.
The mutating variant computes the differential in place of X.
Manifolds.gradient — Methodgradient(M::AbstractManifold, f, p, backend::AbstractRiemannianDiffBackend)
gradient!(M::AbstractManifold, f, X, p, backend::AbstractRiemannianDiffBackend)Compute the Riemannian gradient $∇f(p)$ of a real-valued function $f:\mathcal M \to ℝ$ at point p on the manifold M using the specified AbstractRiemannianDiffBackend.
The mutating variant computes the gradient in place of X.
Manifolds.gradient — Methodgradient(M, f, p, backend::TangentDiffBackend)This method uses the internal backend.diff_backend (Euclidean) on the function
\[ f(\retr_p(\cdot))\]
which is given on the tangent space. In detail, the gradient can be written in terms of the backend.basis_arg. We illustrate it here for an AbstractOrthonormalBasis, since that simplifies notations:
\[\operatorname{grad}f(p) = \operatorname{grad}f(p) = \sum_{i=1}^{d} g_p(\operatorname{grad}f(p),X_i)X_i = \sum_{i=1}^{d} Df(p)[X_i]X_i\]
where the last equality is due to the definition of the gradient as the Riesz representer of the differential.
If the backend is a forward (or backward) finite difference, these coefficients in the sum can be approximates as
\[DF(p)[Y] ≈ \frac{1}{h}\bigl( f(\exp_p(hY)) - f(p) \bigr)\]
writing $p=\exp_p(0)$ we see that this is a finite difference of $f\circ\exp_p$, i.e. for a function on the tangent space, so we can also use other (Euclidean) backends
Manifolds.ExplicitEmbeddedBackend — TypeExplicitEmbeddedBackend{TF<:NamedTuple} <: AbstractDiffBackendA backend to use with the RiemannianProjectionBackend or the TangentDiffBackend, when you have explicit formulae for the gradient in the embedding available.
Constructor
ExplicitEmbeddedBackend(M::AbstractManifold; kwargs)Construct an ExplicitEmbeddedBackend in the embedding M, where currently the following keywords may be used
gradientfor a(n allocating) gradient functiongradient(M, p)defined in the embeddinggradient!for a mutating gradient functiongradient!(M, X, p).
Note that the gradient functions are defined on the embedding manifold M passed to the Backend as well
- AbsilMahonySepulchre2008
Absil, P.-A., Mahony, R. and Sepulchre R., Optimization Algorithms on Matrix Manifolds Princeton University Press, 2008, doi: 10.1515/9781400830244 open access