Generalized Grassmann
Manifolds.GeneralizedGrassmann β TypeGeneralizedGrassmann{n,k,π½} <: AbstractDecoratorManifold{π½}The generalized Grassmann manifold $\operatorname{Gr}(n,k,B)$ consists of all subspaces spanned by $k$ linear independent vectors $π½^n$, where $π½ β \{β, β\}$ is either the real- (or complex-) valued vectors. This yields all $k$-dimensional subspaces of $β^n$ for the real-valued case and all $2k$-dimensional subspaces of $β^n$ for the second.
The manifold can be represented as
\[\operatorname{Gr}(n, k, B) := \bigl\{ \operatorname{span}(p)\ \big|\ p β π½^{n Γ k}, p^\mathrm{H}Bp = I_k\},\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate (or Hermitian) transpose and $I_k$ is the $k Γ k$ identity matrix. This means, that the columns of $p$ form an unitary basis of the subspace with respect to the scaled inner product, that is a point on $\operatorname{Gr}(n,k,B)$, and hence the subspace can actually be represented by a whole equivalence class of representers. For $B=I_n$ this simplifies to the Grassmann manifold.
The tangent space at a point (subspace) $p$ is given by
\[T_x\mathrm{Gr}(n,k,B) = \bigl\{ X β π½^{n Γ k} : X^{\mathrm{H}}Bp + p^{\mathrm{H}}BX = 0_{k} \bigr\},\]
where $0_{k}$ denotes the $k Γ k$ zero matrix.
Note that a point $p β \operatorname{Gr}(n,k,B)$ might be represented by different matrices (i.e. matrices with $B$-unitary column vectors that span the same subspace). Different representations of $p$ also lead to different representation matrices for the tangent space $T_p\mathrm{Gr}(n,k,B)$
The manifold is named after Hermann G. GraΓmann (1809-1877).
Constructor
GeneralizedGrassmann(n, k, B=I_n, field=β)Generate the (real-valued) Generalized Grassmann manifold of $n\times k$ dimensional orthonormal matrices with scalar product B.
Base.exp β Methodexp(M::GeneralizedGrassmann, p, X)Compute the exponential map on the GeneralizedGrassmann M$= \mathrm{Gr}(n,k,B)$ starting in p with tangent vector (direction) X. Let $X^{\mathrm{H}}BX = USV$ denote the SVD decomposition of $X^{\mathrm{H}}BX$. Then the exponential map is written using
\[\exp_p X = p V\cos(S)V^\mathrm{H} + U\sin(S)V^\mathrm{H},\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian and the cosine and sine are applied element wise to the diagonal entries of $S$.
Base.log β Methodlog(M::GeneralizedGrassmann, p, q)Compute the logarithmic map on the GeneralizedGrassmann M$ = \mathcal M=\mathrm{Gr}(n,k,B)$, i.e. the tangent vector X whose corresponding geodesic starting from p reaches q after time 1 on M. The formula reads
\[\log_p q = V\cdot \operatorname{atan}(S) \cdot U^\mathrm{H},\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian. The matrices $U$ and $V$ are the unitary matrices, and $S$ is the diagonal matrix containing the singular values of the SVD-decomposition
\[USV = (q^\mathrm{H}Bp)^{-1} ( q^\mathrm{H} - q^\mathrm{H}Bpp^\mathrm{H}).\]
In this formula the $\operatorname{atan}$ is meant elementwise.
Manifolds.change_metric β Methodchange_metric(M::GeneralizedGrassmann, ::EuclideanMetric, p X)Change X to the corresponding vector with respect to the metric of the GeneralizedGrassmann M, i.e. let $B=LL'$ be the Cholesky decomposition of the matrix M.B, then the corresponding vector is $L\X$.
Manifolds.change_representer β Methodchange_representer(M::GeneralizedGrassmann, ::EuclideanMetric, p, X)Change X to the corresponding representer of a cotangent vector at p with respect to the scaled metric of the GeneralizedGrassmann M, i.e, since
\[g_p(X,Y) = \operatorname{tr}(Y^{\mathrm{H}}BZ) = \operatorname{tr}(X^{\mathrm{H}}Z) = β¨X,Zβ©\]
has to hold for all $Z$, where the repreenter X is given, the resulting representer with respect to the metric on the GeneralizedGrassmann is given by $Y = B^{-1}X$.
ManifoldsBase.check_point β Methodcheck_point(M::GeneralizedGrassmann{n,k,π½}, p)Check whether p is representing a point on the GeneralizedGrassmann M, i.e. its a n-by-k matrix of unitary column vectors with respect to the B inner prudct and of correct eltype with respect to π½.
ManifoldsBase.check_vector β Methodcheck_vector(M::GeneralizedGrassmann{n,k,π½}, p, X; kwargs...)Check whether X is a tangent vector in the tangent space of p on the GeneralizedGrassmann M, i.e. that X is of size and type as well as that
\[ p^{\mathrm{H}}BX + \overline{X^{\mathrm{H}}Bp} = 0_k,\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transpose or Hermitian, $\overline{\cdot}$ the (elementwise) complex conjugate, and $0_k$ denotes the $k Γ k$ zero natrix.
ManifoldsBase.distance β Methoddistance(M::GeneralizedGrassmann, p, q)Compute the Riemannian distance on GeneralizedGrassmann manifold M$= \mathrm{Gr}(n,k,B)$.
Let $USV = p^\mathrm{H}Bq$ denote the SVD decomposition of $p^\mathrm{H}Bq$, where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian. Then the distance is given by
\[d_{\mathrm{Gr}(n,k,B)}(p,q) = \operatorname{norm}(\operatorname{Re}(b)).\]
where
\[b_{i}=\begin{cases} 0 & \text{if} \; S_i β₯ 1\\ \arccos(S_i) & \, \text{if} \; S_i<1. \end{cases}\]
ManifoldsBase.injectivity_radius β Methodinjectivity_radius(M::GeneralizedGrassmann)
injectivity_radius(M::GeneralizedGrassmann, p)Return the injectivity radius on the GeneralizedGrassmann M, which is $\frac{Ο}{2}$.
ManifoldsBase.inner β Methodinner(M::GeneralizedGrassmann, p, X, Y)Compute the inner product for two tangent vectors X, Y from the tangent space of p on the GeneralizedGrassmann manifold M. The formula reads
\[g_p(X,Y) = \operatorname{tr}(X^{\mathrm{H}}BY),\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian.
ManifoldsBase.manifold_dimension β Methodmanifold_dimension(M::GeneralizedGrassmann)Return the dimension of the GeneralizedGrassmann(n,k,π½) manifold M, i.e.
\[\dim \operatorname{Gr}(n,k,B) = k(n-k) \dim_β π½,\]
where $\dim_β π½$ is the real_dimension of π½.
ManifoldsBase.project β Methodproject(M::GeneralizedGrassmann, p, X)Project the n-by-k X onto the tangent space of p on the GeneralizedGrassmann M, which is computed by
\[\operatorname{proj_p}(X) = X - pp^{\mathrm{H}}B^\mathrm{T}X,\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian and $\cdot^{\mathrm{T}}$ the transpose.
ManifoldsBase.project β Methodproject(M::GeneralizedGrassmann, p)Project p from the embedding onto the GeneralizedGrassmann M, i.e. compute q as the polar decomposition of $p$ such that $q^{\mathrm{H}}Bq$ is the identity, where $\cdot^{\mathrm{H}}$ denotes the Hermitian, i.e. complex conjugate transpose.
ManifoldsBase.representation_size β Methodrepresentation_size(M::GeneralizedGrassmann{n,k})Return the represenation size or matrix dimension of a point on the GeneralizedGrassmann M, i.e. $(n,k)$ for both the real-valued and the complex value case.
ManifoldsBase.retract β Methodretract(M::GeneralizedGrassmann, p, X, ::PolarRetraction)Compute the SVD-based retraction PolarRetraction on the GeneralizedGrassmann M, by projecting $p + X$ onto M.
ManifoldsBase.zero_vector β Methodzero_vector(M::GeneralizedGrassmann, p)Return the zero tangent vector from the tangent space at p on the GeneralizedGrassmann M, which is given by a zero matrix the same size as p.
Statistics.mean β Methodmean(
M::GeneralizedGrassmann,
x::AbstractVector,
[w::AbstractWeights,]
method = GeodesicInterpolationWithinRadius(Ο/4);
kwargs...,
)Compute the Riemannian mean of x using GeodesicInterpolationWithinRadius.