Fixed-rank matrices
Manifolds.FixedRankMatrices
โ TypeFixedRankMatrices{m,n,k,๐ฝ} <: Manifold{๐ฝ}
The manifold of $m ร n$ real-valued or complex-valued matrices of fixed rank $k$, i.e.
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 = USV^\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=USV^\mathrm{H}$ is given by
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[Vandereycken2013].
Constructor
FixedRankMatrics(m, n, k[, field=โ])
Generate the manifold of m
-by-n
(field
-valued) matrices of rank k
.
Manifolds.SVDMPoint
โ TypeSVDMPoint <: MPoint
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{T}$ (rows).
Constructors
SVDMPoint(A)
for a matrixA
, stores its svd factors (i.e. implicitly $k=\min\{m,n\}$)SVDMPoint(S)
for anSVD
object, stores its svd factors (i.e. implicitly $k=\min\{m,n\}$)SVDMPoint(U,S,Vt)
for the svd factors to initialize theSVDMPoint
` (i.e. implicitly $k=\min\{m,n\}$)SVDMPoint(A,k)
for a matrixA
, stores its svd factors shortened to the best rank $k$ approximationSVDMPoint(S,k)
for anSVD
object, stores its svd factors shortened to the best rank $k$ approximationSVDMPoint(U,S,Vt,k)
for the svd factors to initialize theSVDMPoint
, stores its svd factors shortened to the best rank $k$ approximation
Manifolds.UMVTVector
โ TypeUMVTVector <: TVector
A tangent vector that can be described as a product $UMV^\mathrm{H}$, at least together with its base point, see for example FixedRankMatrices
. This vector structure stores the additionally (to the point) required fields.
Constructors
UMVTVector(U,M,Vt)
store umv factors to initialize theUMVTVector
UMVTVector(U,M,Vt,k)
store the umv factors after shortening them down to inner dimensions $k$, i.e. in $UMV^\mathrm{H}$, where $M$ is a $k ร k$ matrix.
ManifoldsBase.check_manifold_point
โ Methodcheck_manifold_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
.
ManifoldsBase.check_tangent_vector
โ Methodcheck_tangent_vector(M:FixedRankMatrices{m,n,k}, p, X; check_base_point = true, 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
. The optional parameter check_base_point
indicates, whether to call check_manifold_point
for p
.
ManifoldsBase.inner
โ Methodinner(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
.
ManifoldsBase.manifold_dimension
โ Methodmanifold_dimension(M::FixedRankMatrices{m,n,k,๐ฝ})
Return the manifold dimension for the ๐ฝ
-valued FixedRankMatrices
M
of dimension m
xn
of rank k
, namely
where $\dim_โ ๐ฝ$ is the real_dimension
of ๐ฝ
.
ManifoldsBase.project
โ Methodproject(M, p, A)
project(M, p, X)
Project the matrix $A โ โ^{m,n}$ or a UMVTVector
X
from the embedding or another tangent space onto the tangent space at $p$ on the FixedRankMatrices
M
, further decomposing the result into $X=UMV$, i.e. a UMVTVector
.
ManifoldsBase.representation_size
โ Methodrepresentation_size(M::FixedRankMatrices{m,n,k})
Return the element size of a point on the FixedRankMatrices
M
, i.e. the size of matrices on this manifold $(m,n)$.
ManifoldsBase.retract
โ Methodretract(M, p, X, ::PolarRetraction)
Compute an SVD-based retraction on the FixedRankMatrices
M
by computing
where $U_k S_k V_k^\mathrm{H}$ is the shortened singular value decomposition $USV=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.
ManifoldsBase.zero_tangent_vector
โ Methodzero_tangent_vector(M::FixedRankMatrices, p::SVDMPoint)
Return a UMVTVector
representing the zero tangent vector in the tangent space of p
on the FixedRankMatrices
M
, for example all three elements of the resulting structure are zero matrices.
Literature
- Vandereycken2013
Bart Vandereycken: "Low-rank matrix completion by Riemannian Optimization, SIAM Journal on Optiomoization, 23(2), pp. 1214โ1236, 2013. doi: 10.1137/110845768, arXiv: 1209.3834.