RQ factorization¶
Just like an LQ factorization, but the orthogonalization process starts from the bottom row and produces a much sparser triangular factor when the matrix is wider than it is tall.
Factorization¶
C++ API¶
-
void
RQ
(AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &t, AbstractDistMatrix<Base<F>> &d)¶ Overwrite the matrix \(A\) with \(R\) and the Householder reflectors representing \(\hat Q\). The scalings for the Householder reflectors are stored in the vector t and the unitary diagonal matrix which forces \(R\) to be positive is defined by the vector d.
-
void
rq
::
ExplicitTriang
(AbstractDistMatrix<F> &A)¶ Overwrite A with the triangular factor, R.
C API¶
-
ElError
ElRQDist_z
(ElDistMatrix_z A, ElDistMatrix_z t, ElDistMatrix_d d)¶ Overwrite the matrix \(A\) with \(R\) and the Householder reflectors representing \(\hat Q\). The scalings for the Householder reflectors are stored in the vector t and the unitary diagonal matrix which forces \(R\) to be positive is defined by the vector d.
Apply the factored matrix¶
C++ API¶
-
void
rq
::
ApplyQ
(LeftOrRight side, Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, const Matrix<Base<F>> &d, Matrix<F> &B)¶
-
void
rq
::
ApplyQ
(LeftOrRight side, Orientation orientation, const AbstractDistMatrix<F> &A, const AbstractDistMatrix<F> &t, const AbstractDistMatrix<Base<F>> &d, AbstractDistMatrix<F> &B)¶ Applies the implicitly-defined \(Q\) (or its adjoint) stored within A, t, and d from either the left or the right to \(B\).
C API¶
-
ElError
ElApplyQAfterRQ_s
(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_s A, ElConstMatrix_s t, ElConstMatrix_s d, ElMatrix_s B)¶
-
ElError
ElApplyQAfterRQ_d
(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_d A, ElConstMatrix_d t, ElConstMatrix_d d, ElMatrix_d B)¶
-
ElError
ElApplyQAfterRQ_c
(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_c A, ElConstMatrix_c t, ElConstMatrix_s d, ElMatrix_c B)¶
-
ElError
ElApplyQAfterRQ_z
(ElLeftOrRight side, ElOrientation orientation, ElConstMatrix_z A, ElConstMatrix_z t, ElConstMatrix_d d, ElMatrix_z B)¶
-
ElError
ElApplyQAfterRQDist_s
(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_s A, ElConstDistMatrix_s t, ElConstDistMatrix_s d, ElDistMatrix_s B)¶
-
ElError
ElApplyQAfterRQDist_d
(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_d A, ElConstDistMatrix_d t, ElConstDistMatrix_d d, ElDistMatrix_d B)¶
-
ElError
ElApplyQAfterRQDist_c
(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_c A, ElConstDistMatrix_c t, ElConstDistMatrix_s d, ElDistMatrix_c B)¶
-
ElError
ElApplyQAfterRQDist_z
(ElLeftOrRight side, ElOrientation orientation, ElConstDistMatrix_z A, ElConstDistMatrix_z t, ElConstDistMatrix_d d, ElDistMatrix_z B)¶
Solve linear systems with the factored matrix¶
C++ API¶
-
void
rq
::
SolveAfter
(Orientation orientation, const Matrix<F> &A, const Matrix<F> &t, const Matrix<Base<F>> &d, const Matrix<F> &B, Matrix<F> &X)¶
-
void
rq
::
SolveAfter
(Orientation orientation, const AbstractDistMatrix<F> &A, const AbstractDistMatrix<F> &t, const AbstractDistMatrix<Base<F>> &d, const AbstractDistMatrix<F> &B, AbstractDistMatrix<F> &X)¶ Solves a set of linear systems using an existing packed RQ 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¶
-
ElError
ElSolveAfterRQ_s
(ElOrientation orientation, ElConstMatrix_s A, ElConstMatrix_s t, ElConstMatrix_s d, ElConstMatrix_s B, ElMatrix_s X)¶
-
ElError
ElSolveAfterRQ_d
(ElOrientation orientation, ElConstMatrix_d A, ElConstMatrix_d t, ElConstMatrix_d d, ElConstMatrix_d B, ElMatrix_d X)¶
-
ElError
ElSolveAfterRQ_c
(ElOrientation orientation, ElConstMatrix_c A, ElConstMatrix_c t, ElConstMatrix_s d, ElConstMatrix_c B, ElMatrix_c X)¶
-
ElError
ElSolveAfterRQ_z
(ElOrientation orientation, ElConstMatrix_z A, ElConstMatrix_z t, ElConstMatrix_d d, ElConstMatrix_z B, ElMatrix_z X)¶
-
ElError
ElSolveAfterRQDist_s
(ElOrientation orientation, ElConstDistMatrix_s A, ElConstDistMatrix_s t, ElConstDistMatrix_s d, ElConstDistMatrix_s B, ElDistMatrix_s X)¶
-
ElError
ElSolveAfterRQDist_d
(ElOrientation orientation, ElConstDistMatrix_d A, ElConstDistMatrix_d t, ElConstDistMatrix_d d, ElConstDistMatrix_d B, ElDistMatrix_d X)¶
-
ElError
ElSolveAfterRQDist_c
(ElOrientation orientation, ElConstDistMatrix_c A, ElConstDistMatrix_c t, ElConstDistMatrix_s d, ElConstDistMatrix_c B, ElDistMatrix_c X)¶
-
ElError
ElSolveAfterRQDist_z
(ElOrientation orientation, ElConstDistMatrix_z A, ElConstDistMatrix_z t, ElConstDistMatrix_d d, ElConstDistMatrix_z B, ElDistMatrix_z X)¶