Generalized QR factorization

Implementation

The generalized QR factorization of a pair of matrices \((A,B)\) is analogous to a QR factorization of \(B^{-1} A\) but does not require that \(B\) is square or invertible: unitary matrices \(Q\) and \(Z\), and (right) upper-triangular matrices \(R\) and \(T\), are computed such that

\[A = Q R\]

and

\[B = Q T Z.\]

Thus, if \(B\) was square and invertible, then the QR factorization of \(B^{-1} A\) would be given by \(Z^H (T^{-1} R)\).

C++ API

void GQR(Matrix<F> &A, Matrix<F> &tA, Matrix<Base<F>> &dA, Matrix<F> &B, Matrix<F> &tB, Matrix<Base<F>> &dB)
void GQR(AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &tA, AbstractDistMatrix<Base<F>> &dA, AbstractDistMatrix<F> &B, AbstractDistMatrix<F> &tB, AbstractDistMatrix<Base<F>> &dB)

Overwrite A with both R and the (scaled) Householder vectors which, along with the scalings tA and sign changes dA, define Q. Likewise, B is overwritten with both T and the (scaled) Householder vectors which define Z.

void gqr::ExplicitTriang(Matrix<F> &A, Matrix<F> &B)
void gqr::ExplicitTriang(AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &B)

Overwrite A with R and B with T.

C API

ElError ElGQR_s(ElMatrix_s A, ElMatrix_s tA, ElMatrix_s dA, ElMatrix_s B, ElMatrix_s tB, ElMatrix_s dB)
ElError ElGQR_d(ElMatrix_d A, ElMatrix_d tA, ElMatrix_d dA, ElMatrix_d B, ElMatrix_d tB, ElMatrix_d dB)
ElError ElGQR_c(ElMatrix_c A, ElMatrix_c tA, ElMatrix_s dA, ElMatrix_c B, ElMatrix_c tB, ElMatrix_s dB)
ElError ElGQR_z(ElMatrix_z A, ElMatrix_z tA, ElMatrix_d dA, ElMatrix_z B, ElMatrix_z tB, ElMatrix_d dB)
ElError ElGQRDist_s(ElDistMatrix_s A, ElDistMatrix_s tA, ElDistMatrix_s dA, ElDistMatrix_s B, ElDistMatrix_s tB, ElDistMatrix_s dB)
ElError ElGQRDist_d(ElDistMatrix_d A, ElDistMatrix_d tA, ElDistMatrix_d dA, ElDistMatrix_d B, ElDistMatrix_d tB, ElDistMatrix_d dB)
ElError ElGQRDist_c(ElDistMatrix_c A, ElDistMatrix_c tA, ElDistMatrix_s dA, ElDistMatrix_c B, ElDistMatrix_c tB, ElDistMatrix_s dB)
ElError ElGQRDist_z(ElDistMatrix_z A, ElDistMatrix_z tA, ElDistMatrix_d dA, ElDistMatrix_z B, ElDistMatrix_z tB, ElDistMatrix_d dB)

Overwrite A with both R and the (scaled) Householder vectors which, along with the scalings tA and sign changes dA, define Q. Likewise, B is overwritten with both T and the (scaled) Householder vectors which define Z.

ElError ElGQRExplicitTriang_s(ElMatrix_s A, ElMatrix_s B)
ElError ElGQRExplicitTriang_d(ElMatrix_d A, ElMatrix_d B)
ElError ElGQRExplicitTriang_c(ElMatrix_c A, ElMatrix_c B)
ElError ElGQRExplicitTriang_z(ElMatrix_z A, ElMatrix_z B)
ElError ElGQRExplicitTriangDist_s(ElDistMatrix_s A, ElDistMatrix_s B)
ElError ElGQRExplicitTriangDist_d(ElDistMatrix_d A, ElDistMatrix_d B)
ElError ElGQRExplicitTriangDist_c(ElDistMatrix_c A, ElDistMatrix_c B)
ElError ElGQRExplicitTriangDist_z(ElDistMatrix_z A, ElDistMatrix_z B)

Overwrite A with R and B with T.