Reduction to condensed form

Hermitian to tridiagonal

The currently best-known algorithms for computing eigenpairs of dense Hermitian matrices begin by performing a unitary similarity transformation which reduces the matrix to real symmetric tridiagonal form (usually through Householder transformations). This routine performs said reduction on a Hermitian matrix and stores the scaled Householder vectors in place of the introduced zeroes.

Main source file

Subroutine implementations

Test driver

class HermitianTridiagCtrl
HermitianTridiagApproach approach
GridOrder order
HermitianTridiagCtrl()

Sets approach to HERMITIAN_TRIDIAG_SQUARE and order to ROW_MAJOR.

Note

Please see the Tuning parameters section for extensive information on maximizing the performance of Householder tridiagonalization.

void HermitianTridiag(UpperOrLower uplo, Matrix<F> &A, Matrix<F> &t)
void HermitianTridiag(UpperOrLower uplo, DistMatrix<F> &A, DistMatrix<F, STAR, STAR> &t, const HermitianTridiagCtrl &ctrl = HermitianTridiagCtrl())

Overwrites the main and sub (super) diagonal of the real matrix A with its similar symmetric tridiagonal matrix and stores the scaled Householder vectors below (above) its tridiagonal entries. Complex Hermitian reductions have the added complication of needing to also store the scalings for the Householder vectors (the scaling can be inferred since the Householder vectors must be unit length) if they are to be applied (in the column vector t).

void HermitianTridiag(UpperOrLower uplo, Matrix<F> &A)
void HermitianTridiag(UpperOrLower uplo, DistMatrix<F> &A, const HermitianTridiagCtrl &ctrl = HermitianTridiagCtrl())

Returns just the (appropriate triangle of the) resulting tridiagonal matrix.

hermitian_tridiag namespace

Implementations

void hermitian_tridiag::ApplyQ(LeftOrRight side, UpperOrLower uplo, Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, Matrix<F> &B)
void hermitian_tridiag::ApplyQ(LeftOrRight side, UpperOrLower uplo, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, MD, STAR> &t, DistMatrix<F> &B)
void hermitian_tridiag::ApplyQ(LeftOrRight side, UpperOrLower uplo, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, STAR, STAR> &t, DistMatrix<F> &B)

Apply (from the left or right) the implicitly defined unitary matrix (or its adjoint) represented by the Householder transformations stored within the specified triangle of A and their scalings are stored in the vector t.

Square to Hessenberg

Main header file

Subroutine header files

Test driver

void Hessenberg(UpperOrLower uplo, Matrix<F> &A, Matrix<F> &t)
void Hessenberg(UpperOrLower uplo, DistMatrix<F> &A, DistMatrix<F, STAR, STAR> &t)

Returns the in-place reduction of the matrix A to lower-/upper-Hessenberg form. The vector t contains the scalings for the Householder reflectors, which are stored in the locations of the zeros that they introduced.

void Hessenberg(UpperOrLower uplo, Matrix<F> &A)
void Hessenberg(UpperOrLower uplo, DistMatrix<F> &A)

Returns just the Hessenberg matrix.

hessenberg namespace

TODO

void hessenberg::ApplyQ(UpperOrLower uplo, LeftOrRight side, Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, Matrix<F> &H)
void hessenberg::ApplyQ(UpperOrLower uplo, LeftOrRight side, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, MD, STAR> &t, DistMatrix<F> &H)
void hessenberg::ApplyQ(UpperOrLower uplo, LeftOrRight side, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, STAR, STAR> &t, DistMatrix<F> &H)

General to bidiagonal

Reduces a general fully-populated \(m \times n\) matrix to bidiagonal form through two-sided Householder transformations; when the \(m \ge n\), the result is upper bidiagonal, otherwise it is lower bidiagonal. This routine is most commonly used as a preprocessing step in computing the SVD of a general matrix.

Main header file

Subroutine implementations

Test driver

void Bidiag(Matrix<F> &A, Matrix<F> &tP, Matrix<F> &tQ)
void Bidiag(DistMatrix<F> &A, DistMatrix<F, STAR, STAR> &tP, DistMatrix<F, STAR, STAR> &tQ)

Overwrites the main and sub (or super) diagonal of the real matrix A with the resulting bidiagonal matrix and stores the scaled Householder vectors in the remainder of the matrix. The complex case must also store the scalings of the Householder transformations (in tP and tQ) if they are to be applied.

void Bidiag(Matrix<F> &A)
void Bidiag(DistMatrix<F> &A)

Returns just the resulting bidiagonal matrix.

bidiag namespace

TODO

void bidiag::ApplyQ(LeftOrRight side, Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, Matrix<F> &B)
void bidiag::ApplyQ(LeftOrRight side, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, MD, STAR> &t, DistMatrix<F> &B)
void bidiag::ApplyQ(LeftOrRight side, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, STAR, STAR> &t, DistMatrix<F> &B)
void bidiag::ApplyP(LeftOrRight side, Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, Matrix<F> &B)
void bidiag::ApplyP(LeftOrRight side, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, MD, STAR> &t, DistMatrix<F> &B)
void bidiag::ApplyP(LeftOrRight side, Orientation orientation, const DistMatrix<F> &A, const DistMatrix<F, STAR, STAR> &t, DistMatrix<F> &B)