Generalized Grassmann

Manifolds.GeneralizedGrassmann β€” Type
GeneralizedGrassmann{n,k,𝔽} <: AbstractEmbeddedManifold{𝔽,DefaultEmbeddingType}

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.

source
Base.exp β€” Method
exp(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$.

source
Base.log β€” Method
log(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.

source
ManifoldsBase.check_manifold_point β€” Method
check_manifold_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 𝔽.

source
ManifoldsBase.check_tangent_vector β€” Method
check_tangent_vector(M::GeneralizedGrassmann{n,k,𝔽}, p, X; check_base_point = true, 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. The optional parameter check_base_point indicates, whether to call check_manifold_point for p.

source
ManifoldsBase.distance β€” Method
distance(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}\]
source
ManifoldsBase.inner β€” Method
inner(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.

source
ManifoldsBase.project β€” Method
project(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.

source
ManifoldsBase.project β€” Method
project(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.

source