Fixed-rank matrices

Manifolds.FixedRankMatricesType
FixedRankMatrices{T,𝔽} <: 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 $⋅^{\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.SVDMPointType
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.UMVTVectorType
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.randMethod
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_HessianMethod
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_pointMethod
check_point(M::FixedRankMatrices, p; kwargs...)

Check whether the matrix or SVDMPoint x ids a valid point on the FixedRankMatrices 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_vectorMethod
check_vector(M:FixedRankMatrices, 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.embedMethod
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.embedMethod
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.innerMethod
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.inverse_retractMethod
inverse_retract(M, p, q, ::OrthographicInverseRetraction)

Compute the orthographic inverse retraction FixedRankMatrices M by computing

\[ X = P_{T_{p}M}(q - p) = qVV^\mathrm{T} + UU^{\mathrm{T}}q - UU^{\mathrm{T}}qVV^{\mathrm{T}} - p,\]

where $p$ is a SVDMPoint(U,S,Vt) and $P_{T_{p}M}$ is the projection onto the tangent space at $p$.

For more details, see [AO14].

source
ManifoldsBase.retractMethod
retract(M::FixedRankMatrices, p, X, ::OrthographicRetraction)

Compute the OrthographicRetraction on the FixedRankMatrices M by finding the nearest point to $p + X$ in

\[ p + X + N_{p}\mathcal M \cap \mathcal M\]

where $N_{p}\mathcal M$ is the Normal Space of $T_{p}\mathcal M$.

If $X$ is sufficiently small, then the nearest such point is unique and can be expressed by

\[ q = (U(S + M) + U_{p})(S + M)^{-1}((S + M)V^{\mathrm{T}} + V^{\mathrm{T}}_{p}),\]

where $p$ is a SVDMPoint(U,S,Vt) and $X$ is an UMVTVector(Up,M,Vtp).

For more details, see [AO14].

source
ManifoldsBase.retractMethod
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

[AM12]
P.-A. Absil and J. Malick. Projection-like Retractions on Matrix Manifolds. SIAM Journal on Optimization 22, 135–158 (2012).
[AO14]
P.-A. Absil and I. V. Oseledets. Low-rank retractions: a survey and new results. Computational Optimization and Applications 62, 5–29 (2014).
[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).