Eigensolvers and SVD¶
Hermitian eigensolver¶
Elemental provides a collection of routines for both full and partial solution of the Hermitian eigenvalue problem
where A is the given Hermitian matrix, and unitary Z and real diagonal \(\Lambda\) are sought. In particular, with the eigenvalues and corresponding eigenpairs labeled in non-decreasing order, the three basic modes are:
- Compute all eigenvalues or eigenpairs, \(\{\lambda_i\}_{i=0}^{n-1}\) or \(\{(z_i,\lambda_i)\}_{i=0}^{n-1}\).
- Compute the eigenvalues or eigenpairs with a given range of indices, say \(\{\lambda_i\}_{i=a}^b\) or \(\{(z_i,\lambda_i)\}_{i=a}^b\), with \(0 \le a \le b < n\).
- Compute all eigenpairs (or just eigenvalues) with eigenvalues lying in a particular half-open interval, either \(\{\lambda_i \;|\; \lambda_i \in (a,b] \}\) or \(\{ (z_i,\lambda_i) \;|\; \lambda_i \in (a,b] \}\).
As of now, all three approaches start with Householder tridiagonalization (ala HermitianTridiag()) and then call Matthias Petschow and Paolo Bientinesi’s PMRRR for the tridiagonal eigenvalue problem.
Note
Please see the Tuning parameters section for information on optimizing the reduction to tridiagonal form, as it is the dominant cost in all of Elemental’s Hermitian eigensolvers.
Full spectrum computation¶
- void HermitianEig(UpperOrLower uplo, Matrix<double>& A, Matrix<double>& w)¶
- void HermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<double>& w)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double, VR, STAR>& w)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<double, VR, STAR>& w)¶
Compute the full set of eigenvalues of the double-precision Hermitian matrix A.
- void HermitianEig(UpperOrLower uplo, Matrix<double>& A, Matrix<double>& w, Matrix<double>& Z)¶
- void HermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<double>& w, Matrix<Complex<double>>& Z)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<double, VR, STAR>& w, DistMatrix<Complex<double>>& Z)¶
Compute the full set of eigenpairs of the double-precision Hermitian matrix A.
Index-based subset computation¶
- void HermitianEig(UpperOrLower uplo, Matrix<double>& A, Matrix<double>& w, int a, int b)¶
- void HermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<double>& w, int a, int b)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double, VR, STAR>& w, int a, int b)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<double, VR, STAR>& w, int a, int b)¶
Compute the eigenvalues of a double-precision Hermitian matrix A with indices in the range \(a,a+1,...,b\).
- void HermitianEig(UpperOrLower uplo, Matrix<double>& A, Matrix<double>& w, Matrix<double>& Z, int a, int b)¶
- void HermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<double>& w, Matrix<Complex<double>>& Z)
- void HermitianEig(UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z, int a, int b)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<double, VR, STAR>& w, DistMatrix<Complex<double>>& Z)
Compute the eigenpairs of a double-precision Hermitian matrix A with indices in the range \(a,a+1,...,b\).
Range-based subset computation¶
- void HermitianEig(UpperOrLower uplo, Matrix<double>& A, Matrix<double>& w, double a, double b)¶
- void HermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<double>& w, double a, double b)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double, VR, STAR>& w, double a, double b)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<double, VR, STAR>& w, double a, double b)¶
Compute the eigenvalues of a double-precision Hermitian matrix A lying in the half-open interval \((a,b]\).
- void HermitianEig(UpperOrLower uplo, Matrix<double>& A, Matrix<double>& w, Matrix<double>& Z, double a, double b)¶
- void HermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<double>& w, Matrix<Complex<double>>& Z)
- void HermitianEig(UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z, double a, double b)¶
- void HermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<double, VR, STAR>& w, DistMatrix<Complex<double>>& Z)
Compute the eigenpairs of a double-precision Hermitian matrix A with eigenvalues lying in the half-open interval \((a,b]\).
Sorting the eigenvalues/eigenpairs¶
Since extra time is required in order to sort the eigenvalues/eigenpairs, they are not sorted by default. However, this can be remedied by the appropriate routine from the following list:
- void hermitian_eig::Sort(Matrix<R>& w)¶
- void hermitian_eig::Sort(DistMatrix<R, VR, STAR>& w)¶
Sort a set of eigenvalues in either ascending or descending order.
- void hermitian_eig::Sort(Matrix<typename Base<F>::type>& w, Matrix<F>& Z)¶
- void hermitian_eig::Sort(DistMatrix<typename Base<F>::type, VR, STAR>& w, DistMatrix<F>& Z)¶
Sort a set of eigenpairs in either ascending or descending order (based on the eigenvalues).
Skew-Hermitian eigensolver¶
Essentially identical to the Hermitian eigensolver, HermitianEig(); for any skew-Hermitian matrix \(G\), \(iG\) is Hermitian, as
This fact implies a fast method for solving skew-Hermitian eigenvalue problems:
- Form \(iG\) in \(O(n^2)\) work (switching to complex arithmetic in the real case)
- Run a Hermitian eigensolve on \(iG\), yielding \(iG=Z \Lambda Z^H\).
- Recognize that \(G=Z (-i \Lambda) Z^H\) provides an EVD of \(G\).
Please see the HermitianEig() documentation for more details.
Note
Please see the Tuning parameters section for information on optimizing the reduction to tridiagonal form, as it is the dominant cost in all of Elemental’s Hermitian eigensolvers.
Full spectrum computation¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<double>& G, Matrix<double>& wImag)¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& G, Matrix<double>& wImag)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<double>& G, DistMatrix<double, VR, STAR>& wImag)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& G, DistMatrix<double, VR, STAR>& wImag)¶
Compute the full set of eigenvalues of the double-precision skew-Hermitian matrix G.
- void SkewHermitianEig(UpperOrLower uplo, Matrix<double>& G, Matrix<double>& wImag, Matrix<Complex<double>>& Z)¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& G, Matrix<double>& wImag, Matrix<Complex<double>>& Z)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<double>& G, DistMatrix<double, VR, STAR>& wImag, DistMatrix<Complex<double>>& Z)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& G, DistMatrix<double, VR, STAR>& wImag, DistMatrix<Complex<double>>& Z)¶
Compute the full set of eigenpairs of the double-precision skew-Hermitian matrix G.
Index-based subset computation¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<double>& G, Matrix<double>& wImag, int a, int b)¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& G, Matrix<double>& wImag, int a, int b)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<double>& G, DistMatrix<double, VR, STAR>& wImag, int a, int b)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& G, DistMatrix<double, VR, STAR>& wImag, int a, int b)¶
Compute the eigenvalues of a double-precision skew-Hermitian matrix G with indices in the range \(a,a+1,...,b\).
- void SkewHermitianEig(UpperOrLower uplo, Matrix<double>& G, Matrix<double>& wImag, Matrix<Complex<double>>& Z, int a, int b)¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& G, Matrix<double>& wImag, Matrix<Complex<double>>& Z)
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<double>& G, DistMatrix<double, VR, STAR>& wImag, DistMatrix<Complex<double>>& Z, int a, int b)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& G, DistMatrix<double, VR, STAR>& wImag, DistMatrix<Complex<double>>& Z)
Compute the eigenpairs of a double-precision skew-Hermitian matrix G with indices in the range \(a,a+1,...,b\).
Range-based subset computation¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<double>& G, Matrix<double>& wImag, double a, double b)¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& G, Matrix<double>& wImag, double a, double b)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<double>& G, DistMatrix<double, VR, STAR>& wImag, double a, double b)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& G, DistMatrix<double, VR, STAR>& wImag, double a, double b)¶
Compute the eigenvalues of a double-precision skew-Hermitian matrix G lying in the half-open interval \((a,b]i\).
- void SkewHermitianEig(UpperOrLower uplo, Matrix<double>& G, Matrix<double>& wImag, Matrix<Complex<double>>& Z, double a, double b)¶
- void SkewHermitianEig(UpperOrLower uplo, Matrix<Complex<double>>& G, Matrix<double>& wImag, Matrix<Complex<double>>& Z)
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<double>& G, DistMatrix<double, VR, STAR>& wImag, DistMatrix<Complex<double>>& Z, double a, double b)¶
- void SkewHermitianEig(UpperOrLower uplo, DistMatrix<Complex<double>>& G, DistMatrix<double, VR, STAR>& wImag, DistMatrix<Complex<double>>& Z)
Compute the eigenpairs of a double-precision skew-Hermitian matrix G with eigenvalues lying in the half-open interval \((a,b]i\).
Hermitian generalized-definite eigensolvers¶
The following Hermitian generalized-definite eigenvalue problems frequently appear, where both \(A\) and \(B\) are Hermitian, and \(B\) is additionally positive-definite:
which is denoted with the value ABX via the HermitianGenDefiniteEigType enum,
which uses the BAX value, and finally
which uses the AXBX enum value.
- type HermitianGenDefiniteEigType¶
An enum for specifying either the ABX, BAX, or AXBX generalized eigenvalue problems (described above).
Full spectrum computation¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<double>& A, Matrix<double>& B, Matrix<double>& w)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<Complex<double>>& B, Matrix<double>& w)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double>& B, DistMatrix<double, VR, STAR>& w)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<Complex<double>>& B, DistMatrix<double, VR, STAR>& w)¶
Compute the full set of eigenvalues of a generalized EVP involving the double-precision Hermitian matrices A and B, where B is also positive-definite.
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<double>& A, Matrix<double>& B, Matrix<double>& w, Matrix<double>& Z)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<Complex<double>>& B, Matrix<double>& w, Matrix<double>& Z)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double>& B, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<Complex<double>>& B, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z)¶
Compute the full set of eigenpairs of a generalized EVP involving the double-precision Hermitian matrices A and B, where B is also positive-definite.
Index-based subset computation¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<double>& A, Matrix<double>& B, Matrix<double>& w, int a, int b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<Complex<double>>& B, Matrix<double>& w, int a, int b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double>& B, DistMatrix<double, VR, STAR>& w, int a, int b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<Complex<double>>& B, DistMatrix<double, VR, STAR>& w, int a, int b)¶
Compute the eigenvalues with indices in the range \(a,a+1,...,b\) of a generalized EVP involving the double-precision Hermitian matrices A and B, where B is also positive-definite.
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<double>& A, Matrix<double>& B, Matrix<double>& w, Matrix<double>& Z, int a, int b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<Complex<double>>& B, Matrix<double>& w, Matrix<double>& Z)
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double>& B, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z, int a, int b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<Complex<double>>& B, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z)
Compute the eigenpairs with indices in the range \(a,a+1,...,b\) of a generalized EVP involving the double-precision Hermitian matrices A and B, where B is also positive-definite.
Range-based subset computation¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<double>& A, Matrix<double>& B, Matrix<double>& w, double a, double b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<Complex<double>>& B, Matrix<double>& w, double a, double b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double>& B, DistMatrix<double, VR, STAR>& w, double a, double b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<Complex<double>>& B, DistMatrix<double, VR, STAR>& w, double a, double b)¶
Compute the eigenvalues lying in the half-open interval \((a,b]\) of a generalized EVP involving the double-precision Hermitian matrices A and B, where B is also positive-definite.
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<double>& A, Matrix<double>& B, Matrix<double>& w, Matrix<double>& Z, double a, double b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, Matrix<Complex<double>>& A, Matrix<Complex<double>>& B, Matrix<double>& w, Matrix<double>& Z)
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<double>& A, DistMatrix<double>& B, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z, double a, double b)¶
- void HermitianGenDefiniteEig(HermitianGenDefiniteEigType type, UpperOrLower uplo, DistMatrix<Complex<double>>& A, DistMatrix<Complex<double>>& B, DistMatrix<double, VR, STAR>& w, DistMatrix<double>& Z)
Compute the eigenpairs whose eigenvalues lie in the half-open interval \((a,b]\) of a generalized EVP involving the double-precision Hermitian matrices A and B, where B is also positive-definite.
Unitary eigensolver¶
Not yet written, will likely be based on Ming Gu’s unitary Divide and Conquer algorithm for unitary Hessenberg matrices.
Normal eigensolver¶
Not yet written, will likely be based on Angelika Bunse-Gerstner et al.’s Jacobi-like method for simultaneous diagonalization of the commuting Hermitian and skew-Hermitian portions of the matrix.
Schur decomposition¶
Not yet written, will likely eventually include Greg Henry et al.’s and Robert Granat et al.’s approaches.
Hermitian SVD¶
Given an eigenvalue decomposition of a Hermitian matrix \(A\), say
where \(V\) is unitary and \(\Lambda\) is diagonal and real. Then an SVD of \(A\) can easily be computed as
where the columns of \(U\) equal the columns of \(V\), modulo sign flips introduced by negative eigenvalues.
- void HermitianSVD(UpperOrLower uplo, Matrix<F>& A, Matrix<typename Base<F>::type>& s, Matrix<F>& U, Matrix<F>& V)¶
- void HermitianSVD(UpperOrLower uplo, DistMatrix<F>& A, DistMatrix<typename Base<F>::type, VR, STAR>& s, DistMatrix<F>& U, DistMatrix<F>& V)¶
Return a vector of singular values, \(s\), and the left and right singular vector matrices, \(U\) and \(V\), such that \(A=U \mathrm{diag}(s) V^H\).
- void HermitianSVD(UpperOrLower uplo, Matrix<F>& A, Matrix<typename Base<F>::type>& s)¶
- void HermitianSVD(UpperOrLower uplo, DistMatrix<F>& A, DistMatrix<typename Base<F>::type, VR, STAR>& s)¶
Return the singular values of \(A\) in s. Note that the appropriate triangle of A is overwritten during computation.
Polar decomposition¶
Every matrix \(A\) can be written as \(A=QP\), where \(Q\) is unitary and \(P\) is Hermitian and positive semi-definite. This is known as the polar decomposition of \(A\) and can be constructed as \(Q := U V^H\) and \(P := V \Sigma V^H\), where \(A = U \Sigma V^H\) is the SVD of \(A\). Alternatively, it can be computed through (a dynamically-weighted) Halley iteration.
- void Polar(Matrix<F>& A)¶
- void Polar(DistMatrix<F>& A)¶
- void Polar(Matrix<F>& A, Matrix<F>& P)¶
- void Polar(DistMatrix<F>& A, DistMatrix<F>& P)¶
Compute the polar decomposition of \(A\), \(A=QP\), returning \(Q\) within A and \(P\) within P. The current implementation first computes the SVD.
- void HermitianPolar(UpperOrLower uplo, Matrix<F>& A)¶
- void HermitianPolar(UpperOrLower uplo, DistMatrix<F>& A)¶
- void HermitianPolar(UpperOrLower uplo, Matrix<F>& A, Matrix<F>& P)¶
- void HermitianPolar(UpperOrLower uplo, DistMatrix<F>& A, DistMatrix<F>& P)¶
Compute the polar decomposition through a Hermitian EVD. Since this is equivalent to a Hermitian sign decomposition (if \(\text{sgn}(0)\) is set to 1), these routines are equivalent to HermitianSign.
Detailed interface¶
- int polar::QDWH(Matrix<F>& A, typename Base<F>::type lowerBound, typename Base<F>::type upperBound, int maxits=100 )¶
- int polar::QDWH(DistMatrix<F>& A, typename Base<F>::type lowerBound, typename Base<F>::type upperBound, int maxIts=100 )¶
Overwrites \(A\) with the \(Q\) from the polar decomposition using a QR-based dynamically weighted Halley iteration. The number of iterations used is returned upon completion. TODO: reference to Yuji’s paper
SVD¶
Given a general matrix \(A\), the Singular Value Decomposition is the triplet \((U,\Sigma,V)\) such that
where \(U\) and \(V\) are unitary, and \(\Sigma\) is diagonal with non-negative entries.
- void SVD(Matrix<F>& A, Matrix<typename Base<F>::type>& s, Matrix<F>& V)¶
- void SVD(DistMatrix<F>& A, DistMatrix<typename Base<F>::type, VR, STAR>& s, DistMatrix<F>& V)¶
Overwrites A with \(U\), s with the diagonal entries of \(\Sigma\), and V with \(V\).
- void SVD(Matrix<F>& A, Matrix<typename Base<F>::type>& s)¶
- void SVD(DistMatrix<F>& A, DistMatrix<typename Base<F>::type, VR, STAR>& s)¶
Forms the singular values of \(A\) in s. Note that A is overwritten in order to compute the singular values.