Grassmannian manifold
Manifolds.Grassmann
β TypeGrassmann{n,k,π½} <: AbstractEmbeddedManifold{π½,DefaultIsometricEmbeddingType}
The Grassmann manifold $\operatorname{Gr}(n,k)$ 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) := \bigl\{ \operatorname{span}(p) : p β π½^{n Γ k}, p^\mathrm{H}p = I_k\},\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transpose or Hermitian and $I_k$ is the $k Γ k$ identity matrix. This means, that the columns of $p$ form an unitary basis of the subspace, that is a point on $\operatorname{Gr}(n,k)$, and hence the subspace can actually be represented by a whole equivalence class of representers. Another interpretation is, that
\[\operatorname{Gr}(n,k) = \operatorname{St}(n,k) / \operatorname{O}(k),\]
i.e the Grassmann manifold is the quotient of the Stiefel
manifold and the orthogonal group $\operatorname{O}(k)$ of orthogonal $k Γ k$ matrices.
The tangent space at a point (subspace) $x$ is given by
\[T_x\mathrm{Gr}(n,k) = \bigl\{ X β π½^{n Γ k} : X^{\mathrm{H}}p + p^{\mathrm{H}}X = 0_{k} \bigr\},\]
where $0_k$ is the $k Γ k$ zero matrix.
Note that a point $p β \operatorname{Gr}(n,k)$ might be represented by different matrices (i.e. matrices with 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)$
The manifold is named after Hermann G. GraΓmann (1809-1877).
Constructor
Grassmann(n,k,field=β)
Generate the Grassmann manifold $\operatorname{Gr}(n,k)$, where the real-valued case field = β
is the default.
Base.exp
β Methodexp(M::Grassmann, p, X)
Compute the exponential map on the Grassmann
M
$= \mathrm{Gr}(n,k)$ starting in p
with tangent vector (direction) X
. Let $X = USV$ denote the SVD decomposition of $X$. Then the exponential map is written using
\[z = 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$. A final QR decomposition $z=QR$ is performed for numerical stability reasons, yielding the result as
\[\exp_p X = Q.\]
Base.log
β Methodlog(M::Grassmann, p, q)
Compute the logarithmic map on the Grassmann
M
$ = \mathcal M=\mathrm{Gr}(n,k)$, 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}p)^{-1} ( q^\mathrm{H} - q^\mathrm{H}pp^\mathrm{H}).\]
In this formula the $\operatorname{atan}$ is meant elementwise.
Manifolds.uniform_distribution
β Methoduniform_distribution(M::Grassmann{n,k,β}, p)
Uniform distribution on given (real-valued) Grassmann
M
. Specifically, this is the normalized Haar measure on M
. Generated points will be of similar type as p
.
The implementation is based on Section 2.5.1 in [Chikuse2003]; see also Theorem 2.2.2(iii) in [Chikuse2003].
ManifoldsBase.check_point
β Methodcheck_point(M::Grassmann{n,k,π½}, p)
Check whether p
is representing a point on the Grassmann
M
, i.e. its a n
-by-k
matrix of unitary column vectors and of correct eltype
with respect to π½
.
ManifoldsBase.check_vector
β Methodcheck_vector(M::Grassmann{n,k,π½}, p, X; kwargs...)
Check whether X
is a tangent vector in the tangent space of p
on the Grassmann
M
, i.e. that X
is of size and type as well as that
\[ p^{\mathrm{H}}X + X^{\mathrm{H}}p = 0_k,\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transpose or Hermitian and $0_k$ the $k Γ k$ zero matrix.
ManifoldsBase.distance
β Methoddistance(M::Grassmann, p, q)
Compute the Riemannian distance on Grassmann
manifold M
$= \mathrm{Gr}(n,k)$.
Let $USV = p^\mathrm{H}q$ denote the SVD decomposition of $p^\mathrm{H}q$, where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian. Then the distance is given by
\[d_{\mathrm{Gr}(n,k)}(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::Grassmann)
injectivity_radius(M::Grassmann, p)
Return the injectivity radius on the Grassmann
M
, which is $\frac{Ο}{2}$.
ManifoldsBase.inner
β Methodinner(M::Grassmann, p, X, Y)
Compute the inner product for two tangent vectors X
, Y
from the tangent space of p
on the Grassmann
manifold M
. The formula reads
\[g_p(X,Y) = \operatorname{tr}(X^{\mathrm{H}}Y),\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian.
ManifoldsBase.inverse_retract
β Methodinverse_retract(M::Grassmann, p, q, ::PolarInverseRetraction)
Compute the inverse retraction for the PolarRetraction
, on the Grassmann
manifold M
, i.e.,
\[\operatorname{retr}_p^{-1}q = q*(p^\mathrm{H}q)^{-1} - p,\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian.
ManifoldsBase.inverse_retract
β Methodinverse_retract(M, p, q, ::QRInverseRetraction)
Compute the inverse retraction for the QRRetraction
, on the Grassmann
manifold M
, i.e.,
\[\operatorname{retr}_p^{-1}q = q(p^\mathrm{H}q)^{-1} - p,\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian.
ManifoldsBase.manifold_dimension
β Methodmanifold_dimension(M::Grassmann)
Return the dimension of the Grassmann(n,k,π½)
manifold M
, i.e.
\[\dim \operatorname{Gr}(n,k) = k(n-k) \dim_β π½,\]
where $\dim_β π½$ is the real_dimension
of π½
.
ManifoldsBase.project
β Methodproject(M::Grassmann, p, X)
Project the n
-by-k
X
onto the tangent space of p
on the Grassmann
M
, which is computed by
\[\operatorname{proj_p}(X) = X - pp^{\mathrm{H}}X,\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian.
ManifoldsBase.representation_size
β Methodrepresentation_size(M::Grassmann{n,k})
Return the represenation size or matrix dimension of a point on the Grassmann
M
, i.e. $(n,k)$ for both the real-valued and the complex value case.
ManifoldsBase.retract
β Methodretract(M::Grassmann, p, X, ::PolarRetraction)
Compute the SVD-based retraction PolarRetraction
on the Grassmann
M
. With $USV = p + X$ the retraction reads
\[\operatorname{retr}_p X = UV^\mathrm{H},\]
where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transposed or Hermitian.
ManifoldsBase.retract
β Methodretract(M::Grassmann, p, X, ::QRRetraction )
Compute the QR-based retraction QRRetraction
on the Grassmann
M
. With $QR = p + X$ the retraction reads
\[\operatorname{retr}_p X = QD,\]
where D is a $m Γ n$ matrix with
\[D = \operatorname{diag}\left( \operatorname{sgn}\left(R_{ii}+\frac{1}{2}\right)_{i=1}^n \right).\]
ManifoldsBase.vector_transport_to
β Methodvector_transport_to(M::Grassmann,p,X,q,::ProjectionTransport)
compute the projection based transport on the Grassmann
M
by interpreting X
from the tangent space at p
as a point in the embedding and projecting it onto the tangent space at q.
ManifoldsBase.zero_vector
β Methodzero_vector(M::Grassmann, p)
Return the zero tangent vector from the tangent space at p
on the Grassmann
M
, which is given by a zero matrix the same size as p
.
Statistics.mean
β Methodmean(
M::Grassmann,
x::AbstractVector,
[w::AbstractWeights,]
method = GeodesicInterpolationWithinRadius(Ο/4);
kwargs...,
)
Compute the Riemannian mean
of x
using GeodesicInterpolationWithinRadius
.
- Chikuse2003
Y. Chikuse: "Statistics on Special Manifolds", Springer New York, 2003, doi: 10.1007/978-0-387-21540-2.