Convex optimization¶
LogBarrier¶
Uses a careful calculation of the log of the determinant in order to return the log barrier of a Hermitian positive-definite matrix A, \(-\log(\mbox{det}(A))\).
-
Base<F>
LogBarrier
(UpperOrLower uplo, const Matrix<F> &A)¶
-
Base<F>
LogBarrier
(UpperOrLower uplo, const DistMatrix<F> &A)¶
-
Base<F>
LogBarrier
(UpperOrLower uplo, Matrix<F> &A, bool canOverwrite = false)¶
-
Base<F>
LogBarrier
(UpperOrLower uplo, DistMatrix<F> &A, bool canOverwrite = false)¶
LogDetDivergence¶
The log-det divergence of a pair of \(n \times n\) Hermitian positive-definite matrices \(A\) and \(B\) is
which can be greatly simplified using the Cholesky factors of \(A\) and \(B\). In particular, if we set \(Z = L_B^{-1} L_A\), where \(A=L_A L_A^H\) and \(B=L_B L_B^H\) are Cholesky factorizations, then
-
Base<F>
LogDetDivergence
(UpperOrLower uplo, const Matrix<F> &A, const Matrix<F> &B)¶
-
Base<F>
LogDetDivergence
(UpperOrLower uplo, const DistMatrix<F> &A, const DistMatrix<F> &B)¶
Singular-value soft-thresholding¶
Overwrites \(A\) with \(U S_{\tau}(\Sigma) V^H\), where \(U \Sigma V^H\) is the singular-value decomposition of \(A\) upon input and \(S_{\tau}\) performs soft-thresholding with parameter \(\tau\). The return value is the rank of the soft-thresholded matrix.
-
int
SVT
(DistMatrix<F> &A, Base<F> tau, bool relative = false)¶ Runs the default SVT algorithm. In the sequential case, this is currently svt::Normal, and, in the parallel case, it is svt::Cross.
-
int
SVT
(DistMatrix<F> &A, Base<F> tau, int relaxedRank, bool relative = false)¶ Runs a faster (for small ranks), but less accurate, algorithm given an upper bound on the rank of the soft-thresholded matrix. The current implementation preprocesses via relaxedRank steps of (Businger-Golub) column-pivoted QR via the routine svt::PivotedQR.
-
int
SVT
(DistMatrix<F, U, STAR> &A, Base<F> tau, bool relative = false)¶ Runs an SVT algorithm designed for tall-skinny matrices. The current implementation is based on TSQR factorization and is svt::TSQR.
namespace svt¶
-
int
svt
::
Normal
(DistMatrix<F> &A, Base<F> tau, bool relative = false)¶ Runs a standard SVD, soft-thresholds the singular values, and then reforms the matrix.
-
int
svt
::
Cross
(DistMatrix<F> &A, Base<F> tau, bool relative = false)¶ Forms the normal matrix, computes its Hermitian EVD, soft-thresholds the eigenvalues, and then reforms the matrix. Note that Elemental’s parallel Hermitian EVD is much faster than its parallel SVD; this is typically worth the loss of accuracy in the computed small (truncated) singular values and is therefore the default choice for parallel SVT.
-
int
svt
::
PivotedQR
(DistMatrix<F> &A, Base<F> tau, int numStepsQR, bool relative = false)¶ Computes an approximate SVT by first approximating A as the rank-numSteps approximation produced by numSteps iterations of column-pivoted QR.