LQ factorization

Implementation

Subroutines

Test driver

Given \(A \in \mathbb{F}^{m \times n}\), an LQ factorization typically computes an implicit unitary matrix \(\hat Q \in \mathbb{F}^{n \times n}\) such that \(\hat L \equiv A\hat Q^H\) is lower trapezoidal. One can then form the thin factors \(L \in \mathbb{F}^{m \times \mbox{min}(m,n)}\) and \(Q \in \mathbb{F}^{\mbox{min}(m,n) \times n}\) by setting \(L\) and \(Q\) to first \(\mbox{min}(m,n)\) columns and rows of \(\hat L\) and \(\hat Q\), respectively. Upon completion \(L\) is stored in the lower trapezoid of \(A\) and the Householder reflectors (and preceding unitary diagonal matrix forcing \(L\) to have a positive diagonal, defined by the vector d) representing \(\hat Q\) are stored within the rows of the strictly upper trapezoid.

Factorization

C++ API

void LQ(Matrix<F> &A, Matrix<F> &t, Matrix<Base<F>> &d)
void LQ(AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &t, AbstractDistMatrix<Base<F>> &d)

Overwrite the matrix \(A\) with \(L\) and the Householder reflectors representing \(\hat Q\). The scalings for the Householder reflectors are stored in the vector t and the diagonal matrix which forces \(L\) to be positive in d.

void lq::Explicit(Matrix<F> &L, Matrix<F> &A)
void lq::Explicit(AbstractDistMatrix<F> &L, AbstractDistMatrix<F> &A)

Overwrite A with Q and return L.

void lq::ExplicitTriang(Matrix<F> &A)
void lq::ExplicitTriang(AbstractDistMatrix<F> &A)

Overwrite A with the triangular factor, L.

void lq::ExplicitUnitary(Matrix<F> &A)
void lq::ExplicitUnitary(AbstractDistMatrix<F> &A)

Overwrite A with Q.

C API

ElError ElLQ_s(ElMatrix_s A, ElMatrix_s t, ElMatrix_s d)
ElError ElLQ_d(ElMatrix_d A, ElMatrix_d t, ElMatrix_d d)
ElError ElLQ_c(ElMatrix_c A, ElMatrix_c t, ElMatrix_s d)
ElError ElLQ_z(ElMatrix_z A, ElMatrix_z t, ElMatrix_d d)
ElError ElLQDist_s(ElDistMatrix_s A, ElDistMatrix_s t, ElDistMatrix_s d)
ElError ElLQDist_d(ElDistMatrix_d A, ElDistMatrix_d t, ElDistMatrix_d d)
ElError ElLQDist_c(ElDistMatrix_c A, ElDistMatrix_c t, ElDistMatrix_s d)
ElError ElLQDist_z(ElDistMatrix_z A, ElDistMatrix_z t, ElDistMatrix_d d)

Overwrite the matrix \(A\) with \(L\) and the Householder reflectors representing \(\hat Q\). The scalings for the Householder reflectors are stored in the vector t and the diagonal matrix which forces \(L\) to be positive in d.

ElError ElLQExplicit_s(ElMatrix_s L, ElMatrix_s A)
ElError ElLQExplicit_d(ElMatrix_d L, ElMatrix_d A)
ElError ElLQExplicit_c(ElMatrix_c L, ElMatrix_c A)
ElError ElLQExplicit_z(ElMatrix_z L, ElMatrix_z A)
ElError ElLQExplicitDist_s(ElDistMatrix_s L, ElDistMatrix_s A)
ElError ElLQExplicitDist_d(ElDistMatrix_d L, ElDistMatrix_d A)
ElError ElLQExplicitDist_c(ElDistMatrix_c L, ElDistMatrix_c A)
ElError ElLQExplicitDist_z(ElDistMatrix_z L, ElDistMatrix_z A)

Overwrite A with Q and return L.

ElError ElLQExplicitTriang_s(ElMatrix_s A)
ElError ElLQExplicitTriang_d(ElMatrix_d A)
ElError ElLQExplicitTriang_c(ElMatrix_c A)
ElError ElLQExplicitTriang_z(ElMatrix_z A)
ElError ElLQExplicitTriangDist_s(ElDistMatrix_s A)
ElError ElLQExplicitTriangDist_d(ElDistMatrix_d A)
ElError ElLQExplicitTriangDist_c(ElDistMatrix_c A)
ElError ElLQExplicitTriangDist_z(ElDistMatrix_z A)

Ovewrite A with the triangular factor, L.

ElError ElLQExplicitUnitary_s(ElMatrix_s A)
ElError ElLQExplicitUnitary_d(ElMatrix_d A)
ElError ElLQExplicitUnitary_c(ElMatrix_c A)
ElError ElLQExplicitUnitary_z(ElMatrix_z A)
ElError ElLQExplicitUnitaryDist_s(ElDistMatrix_s A)
ElError ElLQExplicitUnitaryDist_d(ElDistMatrix_d A)
ElError ElLQExplicitUnitaryDist_c(ElDistMatrix_c A)
ElError ElLQExplicitUnitaryDist_z(ElDistMatrix_z A)

Overwrite A with Q.

Applying the factored matrix

The following routines apply the implicitly-defined \(Q\) (or its adjoint) stored within A, t, and d from either the left or the right to \(B\).

C++ API

void lq::ApplyQ(LeftOrRight side, Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, const Matrix<Base<F>> &d, Matrix<F> &B)
void lq::ApplyQ(LeftOrRight side, Orientation orientation, const AbstractDistMatrix<F> &A, const AbstractDistMatrix<F> &t, const AbstractDistMatrix<Base<F>> &d, AbstractDistMatrix<F> &B)

C API

ElError ElApplyQAfterLQ_s(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_s A, ElConstMatrix_s t, ElConstMatrix_s d, ElMatrix_s B)
ElError ElApplyQAfterLQ_d(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_d A, ElConstMatrix_d t, ElConstMatrix_d d, ElMatrix_d B)
ElError ElApplyQAfterLQ_c(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_c A, ElConstMatrix_c t, ElConstMatrix_s d, ElMatrix_c B)
ElError ElApplyQAfterLQ_z(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_z A, ElConstMatrix_z t, ElConstMatrix_d d, ElMatrix_z B)
ElError ElApplyQAfterLQDist_s(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_s A, ElConstDistMatrix_s t, ElConstDistMatrix_s d, ElDistMatrix_s B)
ElError ElApplyQAfterLQDist_d(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_d A, ElConstDistMatrix_d t, ElConstDistMatrix_d d, ElDistMatrix_d B)
ElError ElApplyQAfterLQDist_c(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_c A, ElConstDistMatrix_c t, ElConstDistMatrix_s d, ElDistMatrix_c B)
ElError ElApplyQAfterLQDist_z(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_z A, ElConstDistMatrix_z t, ElConstDistMatrix_d d, ElDistMatrix_z B)

Solving against the factored matrix

The following routines solve a set of linear systems using an existing packed LQ factorization given by \(A\) and the vectors \(t\) and \(d\). \(B\) is the matrix of input vectors and \(X\) is the matrix of solutions.

C++ API

void lq::SolveAfter(Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, const Matrix<Base<F>> &d, const Matrix<F> &B, Matrix<F> &X)
void lq::SolveAfter(Orientation orientation, const AbstractDistMatrix<F> &A, const AbstractDistMatrix<F> &t, const AbstractDistMatrix<Base<F>> &d, const AbstractDistMatrix<F> &B, AbstractDistMatrix<F> &X)

C API

ElError ElSolveAfterLQ_s(ElOrientation orientation, ElConstMatrix_s A, ElConstMatrix_s t, ElConstMatrix_s d, ElConstMatrix_s B, ElMatrix_s X)
ElError ElSolveAfterLQ_d(ElOrientation orientation, ElConstMatrix_d A, ElConstMatrix_d t, ElConstMatrix_d d, ElConstMatrix_d B, ElMatrix_d X)
ElError ElSolveAfterLQ_c(ElOrientation orientation, ElConstMatrix_c A, ElConstMatrix_c t, ElConstMatrix_s d, ElConstMatrix_c B, ElMatrix_c X)
ElError ElSolveAfterLQ_z(ElOrientation orientation, ElConstMatrix_z A, ElConstMatrix_z t, ElConstMatrix_d d, ElConstMatrix_z B, ElMatrix_z X)
ElError ElSolveAfterLQDist_s(ElOrientation orientation, ElConstDistMatrix_s A, ElConstDistMatrix_s t, ElConstDistMatrix_s d, ElConstDistMatrix_s B, ElDistMatrix_s X)
ElError ElSolveAfterLQDist_d(ElOrientation orientation, ElConstDistMatrix_d A, ElConstDistMatrix_d t, ElConstDistMatrix_d d, ElConstDistMatrix_d B, ElDistMatrix_d X)
ElError ElSolveAfterLQDist_c(ElOrientation orientation, ElConstDistMatrix_c A, ElConstDistMatrix_c t, ElConstDistMatrix_s d, ElConstDistMatrix_c B, ElDistMatrix_c X)
ElError ElSolveAfterLQDist_z(ElOrientation orientation, ElConstDistMatrix_z A, ElConstDistMatrix_z t, ElConstDistMatrix_d d, ElConstDistMatrix_z B, ElDistMatrix_z X)