Fixed-rank matrices

Manifolds.FixedRankMatrices โ€” Type
FixedRankMatrices{m,n,k,๐”ฝ} <: AbstractDecoratorManifold{๐”ฝ}

The manifold of $m ร— n$ real-valued or complex-valued matrices of fixed rank $k$, i.e.

\[\bigl\{ p โˆˆ ๐”ฝ^{m ร— n}\ \big|\ \operatorname{rank}(p) = k\bigr\},\]

where $๐”ฝ โˆˆ \{โ„,โ„‚\}$ and the rank is the number of linearly independent columns of a matrix.

Representation with 3 matrix factors

A point $p โˆˆ \mathcal M$ can be stored using unitary matrices $U โˆˆ ๐”ฝ^{m ร— k}$, $V โˆˆ ๐”ฝ^{n ร— k}$ as well as the $k$ singular values of $p = U_p S V_p^\mathrm{H}$, where $\cdot^{\mathrm{H}}$ denotes the complex conjugate transpose or Hermitian. In other words, $U$ and $V$ are from the manifolds Stiefel(m,k,๐”ฝ) and Stiefel(n,k,๐”ฝ), respectively; see SVDMPoint for details.

The tangent space $T_p \mathcal M$ at a point $p โˆˆ \mathcal M$ with $p=U_p S V_p^\mathrm{H}$ is given by

\[T_p\mathcal M = \bigl\{ U_p M V_p^\mathrm{H} + U_X V_p^\mathrm{H} + U_p V_X^\mathrm{H} : M โˆˆ ๐”ฝ^{k ร— k}, U_X โˆˆ ๐”ฝ^{m ร— k}, V_X โˆˆ ๐”ฝ^{n ร— k} \text{ s.t. } U_p^\mathrm{H}U_X = 0_k, V_p^\mathrm{H}V_X = 0_k \bigr\},\]

where $0_k$ is the $k ร— k$ zero matrix. See UMVTVector for details.

The (default) metric of this manifold is obtained by restricting the metric on $โ„^{m ร— n}$ to the tangent bundle [Van13].

Constructor

FixedRankMatrices(m, n, k[, field=โ„])

Generate the manifold of m-by-n (field-valued) matrices of rank k.

source
Manifolds.SVDMPoint โ€” Type
SVDMPoint <: AbstractManifoldPoint

A point on a certain manifold, where the data is stored in a svd like fashion, i.e. in the form $USV^\mathrm{H}$, where this structure stores $U$, $S$ and $V^\mathrm{H}$. The storage might also be shortened to just $k$ singular values and accordingly shortened $U$ (columns) and $V^\mathrm{H}$ (rows).

Constructors

  • SVDMPoint(A) for a matrix A, stores its svd factors (i.e. implicitly $k=\min\{m,n\}$)
  • SVDMPoint(S) for an SVD object, stores its svd factors (i.e. implicitly $k=\min\{m,n\}$)
  • SVDMPoint(U,S,Vt) for the svd factors to initialize the SVDMPoint(i.e. implicitlyk=\min\{m,n\}`)
  • SVDMPoint(A,k) for a matrix A, stores its svd factors shortened to the best rank $k$ approximation
  • SVDMPoint(S,k) for an SVD object, stores its svd factors shortened to the best rank $k$ approximation
  • SVDMPoint(U,S,Vt,k) for the svd factors to initialize the SVDMPoint, stores its svd factors shortened to the best rank $k$ approximation
source
Manifolds.UMVTVector โ€” Type
UMVTVector <: TVector

A tangent vector that can be described as a product $U_p M V_p^\mathrm{H} + U_X V_p^\mathrm{H} + U_p V_X^\mathrm{H}$, where $X = U_X S V_X^\mathrm{H}$ is its base point, see for example FixedRankMatrices.

The base point $p$ is required for example embedding this point, but it is not stored. The fields of thie tangent vector are U for $U_X$, M and Vt to store $V_X^\mathrm{H}$

Constructors

  • UMVTVector(U,M,Vt) store umv factors to initialize the UMVTVector
  • UMVTVector(U,M,Vt,k) store the umv factors after shortening them down to inner dimensions k.
source
Base.rand โ€” Method
Random.rand(M::FixedRankMatrices; vector_at=nothing, kwargs...)

If vector_at is nothing, return a random point on the FixedRankMatrices manifold. The orthogonal matrices are sampled from the Stiefel manifold and the singular values are sampled uniformly at random.

If vector_at is not nothing, generate a random tangent vector in the tangent space of the point vector_at on the FixedRankMatrices manifold M.

source
ManifoldDiff.riemannian_Hessian โ€” Method
Y = riemannian_Hessian(M::FixedRankMatrices, p, G, H, X)
riemannian_Hessian!(M::FixedRankMatrices, Y, p, G, H, X)

Compute the Riemannian Hessian $\operatorname{Hess} f(p)[X]$ given the Euclidean gradient $โˆ‡ f(\tilde p)$ in G and the Euclidean Hessian $โˆ‡^2 f(\tilde p)[\tilde X]$ in H, where $\tilde p, \tilde X$ are the representations of $p,X$ in the embedding,.

The Riemannian Hessian can be computed as stated in Remark 4.1 [Ngu23] or Section 2.3 [Van13], that B. Vandereycken adopted for Manopt (Matlab).

source
ManifoldsBase.check_point โ€” Method
check_point(M::FixedRankMatrices{m,n,k}, p; kwargs...)

Check whether the matrix or SVDMPoint x ids a valid point on the FixedRankMatrices{m,n,k,๐”ฝ} M, i.e. is an m-byn matrix of rank k. For the SVDMPoint the internal representation also has to have the right shape, i.e. p.U and p.Vt have to be unitary. The keyword arguments are passed to the rank function that verifies the rank of p.

source
ManifoldsBase.check_vector โ€” Method
check_vector(M:FixedRankMatrices{m,n,k}, p, X; kwargs...)

Check whether the tangent UMVTVector X is from the tangent space of the SVDMPoint p on the FixedRankMatrices M, i.e. that v.U and v.Vt are (columnwise) orthogonal to x.U and x.Vt, respectively, and its dimensions are consistent with p and X.M, i.e. correspond to m-by-n matrices of rank k.

source
ManifoldsBase.embed โ€” Method
embed(M::FixedRankMatrices, p, X)

Embed the tangent vector X at point p in M from its UMVTVector representation into the set of $mร—n$ matrices.

The formula reads

\[U_pMV_p^{\mathrm{H}} + U_XV_p^{\mathrm{H}} + U_pV_X^{\mathrm{H}}\]

source
ManifoldsBase.embed โ€” Method
embed(::FixedRankMatrices, p::SVDMPoint)

Embed the point p from its SVDMPoint representation into the set of $mร—n$ matrices by computing $USV^{\mathrm{H}}$.

source
ManifoldsBase.inner โ€” Method
inner(M::FixedRankMatrices, p::SVDMPoint, X::UMVTVector, Y::UMVTVector)

Compute the inner product of X and Y in the tangent space of p on the FixedRankMatrices M, which is inherited from the embedding, i.e. can be computed using dot on the elements (U, Vt, M) of X and Y.

source
ManifoldsBase.manifold_dimension โ€” Method
manifold_dimension(M::FixedRankMatrices{m,n,k,๐”ฝ})

Return the manifold dimension for the ๐”ฝ-valued FixedRankMatrices M of dimension mxn of rank k, namely

\[\dim(\mathcal M) = k(m + n - k) \dim_โ„ ๐”ฝ,\]

where $\dim_โ„ ๐”ฝ$ is the real_dimension of ๐”ฝ.

source
ManifoldsBase.retract โ€” Method
retract(M, p, X, ::PolarRetraction)

Compute an SVD-based retraction on the FixedRankMatrices M by computing

\[ q = U_kS_kV_k^\mathrm{H},\]

where $U_k S_k V_k^\mathrm{H}$ is the shortened singular value decomposition $USV^\mathrm{H}=p+X$, in the sense that $S_k$ is the diagonal matrix of size $k ร— k$ with the $k$ largest singular values and $U$ and $V$ are shortened accordingly.

source

Literature

[HU17]
S. Hosseini and A. Uschmajew. A Riemannian Gradient Sampling Algorithm for Nonsmooth Optimization on Manifolds. SIAM J. Optim. 27, 173โ€“189 (2017).
[Ngu23]
D. Nguyen. Operator-Valued Formulas for Riemannian Gradient and Hessian and Families of Tractable Metrics in Riemannian Optimization. Journal of Optimization Theory and Applications 198, 135โ€“164 (2023), arXiv:2009.10159.
[Van13]
B. Vandereycken. Low-rank matrix completion by Riemannian optimization. SIAM Journal on Optimization 23, 1214โ€“1236 (2013).