Orthogonal Factorizations

Driver Routines

template<Number T, DataLocation L, template<typename, typename> typename S>
QR<T, L, kInPlace> qr(const Matrix<T, L, S> &A, InPlace)

Factorizes the matrix A as \( \mathbf{A} = \mathbf{QR}\), where \( \mathbf{Q} \) is a unitary matrix and \( \mathbf{R} \) is an upper triangular matrix.

Warning

The factorization is done in-place, overwriting the matrix A.

Parameters:

A[in] input matrix

Returns:

An object representing the factorization.

template<Number T, DataLocation L, template<typename, typename> typename S>
QR<T, L> qr(const Matrix<T, L, S> &A)

Factorizes the matrix A as \( \mathbf{A} = \mathbf{QR}\), where \( \mathbf{Q} \) is a unitary matrix and \( \mathbf{R} \) is an upper triangular matrix.

Parameters:

A[in] input matrix

Returns:

An object representing the factorization.

template<Number T, DataLocation L, template<typename, typename> typename S>
QR<T, L, kInPlace | kColPivot> pivoted_qr(const Matrix<T, L, S> &A, InPlace)

Factorizes the matrix A as \( \mathbf{AP} = \mathbf{QR}\), where \( \mathbf{Q} \) is a unitary matrix, \( \mathbf{R} \) is an upper triangular matrix and \( \mathbf{P} \) is a permutation matrix.

Warning

The factorization is done in-place, overwriting the matrix A.

Parameters:

A[in] input matrix

Returns:

An object representing the factorization.

template<Number T, DataLocation L, template<typename, typename> typename S>
QR<T, L, kColPivot> pivoted_qr(const Matrix<T, L, S> &A)

Factorizes the matrix A as \( \mathbf{AP} = \mathbf{QR}\), where \( \mathbf{Q} \) is a unitary matrix, \( \mathbf{R} \) is an upper triangular matrix and \( \mathbf{P} \) is a permutation matrix.

Parameters:

A[in] input matrix

Returns:

An object representing the factorization.

template<Number T, DataLocation L, template<typename, typename> typename S>
LQ<T, L, kInPlace> lq(const Matrix<T, L, S> &A, InPlace)

Factorizes the matrix A as \( \mathbf{A} = \mathbf{LQ}\), where \( \mathbf{Q} \) is a unitary matrix and \( \mathbf{L} \) is a lower triangular matrix.

Warning

The factorization is done in-place, overwriting the matrix A.

Parameters:

A[in] input matrix

Returns:

An object representing the factorization.

template<Number T, DataLocation L, template<typename, typename> typename S>
LQ<T, L> lq(const Matrix<T, L, S> &A)

Factorizes the matrix A as \( \mathbf{A} = \mathbf{LQ}\), where \( \mathbf{Q} \) is a unitary matrix and \( \mathbf{R} \) is a lower triangular matrix.

Parameters:

A[in] input matrix

Returns:

An object representing the factorization.

template<typename T, typename L, template<typename, typename> typename S1, template<typename, typename> typename S2>
internal::QRSolver<T, L, Matrix<T, L, S2>> qrsolve(const Matrix<T, L, S1> &A, const Matrix<T, L, S2> &B)

Solves the underdetermined or overdetermined linear system \( \mathbf{A}\mathbf{X} = \mathbf{B}\) using QR or LQ factorization such that \( \min_{\mathbf{X}} \| \mathbf{A}\mathbf{X} - \mathbf{B} \|\) and the matrix A has full rank.

Parameters:
  • A[in] coefficient matrix

  • B[in] right-hand side of the system

Returns:

An abstract object representing the solution of the linear system.

template<typename T, typename L, template<typename, typename> typename S1, template<typename, typename> typename S2>
internal::QRSolver<T, L, Matrix<T, L, S2>, kInPlace> qrsolve(const Matrix<T, L, S1> &A, const Matrix<T, L, S2> &B, InPlace)

Solves the underdetermined or overdetermined linear system \( \mathbf{A}\mathbf{X} = \mathbf{B}\) using QR or LQ factorization such that \( \min_{\mathbf{X}} \| \mathbf{A}\mathbf{X} - \mathbf{B} \|\) and the matrix A has full rank.

Warning

The matrix A will be overwritten by the factor Q and R or L of the factorization.

Parameters:
  • A[in] coefficient matrix

  • B[in] right-hand side of the system

Returns:

An abstract object representing the solution of the linear system.

template<typename T, typename L, template<typename, typename> typename S1, template<typename, typename> typename S2>
internal::QRSolver<T, L, Vector<T, L, S2>> qrsolve(const Matrix<T, L, S1> &A, const Vector<T, L, S2> &b)

Solves the underdetermined or overdetermined linear system \( \mathbf{A}\vec{x} = \vec{b}\) using QR or LQ factorization such that \( \min_{\vec{x}} \| \mathbf{A}\vec{x} - \vec{b} \|\) and the matrix A has full rank.

Parameters:
  • A[in] coefficient matrix

  • b[in] right-hand side of the system

Returns:

An abstract object representing the solution of the linear system.

template<typename T, typename L, template<typename, typename> typename S1, template<typename, typename> typename S2>
internal::QRSolver<T, L, Vector<T, L, S2>> qrsolve(const Matrix<T, L, S1> &A, const Vector<T, L, S2> &b, InPlace)

Solves the underdetermined or overdetermined linear system \( \mathbf{A}\vec{x} = \vec{b}\) using QR or LQ factorization such that \( \min_{\vec{x}} \| \mathbf{A}\vec{x} - \vec{b} \|\) and the matrix A has full rank.

Warning

The matrix A will be overwritten by the factor Q and R or L of the factorization.

Parameters:
  • A[in] coefficient matrix

  • b[in] right-hand side of the system

Returns:

An abstract object representing the solution of the linear system.

QR

template<Number T, DataLocation L, int Flags>
class QR

A solver representing the QR factorization of an \( m \times n \) matrix A:

\[ \mathbf{A} = \mathbf{QR} \]
where Q is an orthogonal matrix and R is an upper triangular.

If \( m \geq n \), then the QR factorization is given by

\[\begin{split} \mathbf{A} = \begin{bmatrix} \mathbf{Q_1} & \mathbf{Q_2} \end{bmatrix} \begin{bmatrix} \mathbf{R_1} \\ 0 \end{bmatrix} \end{split}\]
where R is an \( n \times n \) upper triangular matrix with real diagonal elements, and Q is an \( m \times m \) orthogonal (or unitary) matrix. The matrix \( \mathbf{Q_1} \mathbf{R_1} \) is referred as the thin QR factorization of the matrix A.

Usage:

Matrix<double> A(n, n);
Matrix<double> B(n, n);
Matrix<double> Q(n, n);
QR<double> F = qr(A);
Q = F.Q();
B = A * F.Q();

Note

Multiplying the matrix qr.thin_Q() directly as shown by the example is only allowed when m <= n. Otherwise, you need to form the matrix Q1 explicitly: Q1 = qr.thin_Q().

Note

When using column-pivoting (Flags & colPivot), this solver computes \( \mathbf{AP} = \mathbf{QR} \), where \( \mathbf{P} \) is a permutation matrix. Use the reorder method to return to the original problem.

Template Parameters:
  • T – Floating-point type

  • L – Location of the data.

  • Flags – Flags for the QR decomposition. Available: kDefaultFlags (default flags), kInPlace (the factorization is done in-place, overwriting A) and colPivot (enables column pivoting).

Public Functions

inline QR()

Default constructor.

template<template<typename, typename> typename S>
inline QR(const Matrix<T, L, S> &A)

Creates a new solver and factorize the matrix A.

Parameters:

A[in] input matrix

virtual ~QR() = default

Destructor.

template<template<typename, typename> typename S>
inline QR &factorize(const Matrix<T, L, S> &A)

Calculates the QR factorization of the matrix A.

Parameters:

A[in] input matrix

inline HouseholderQ<T, L> Q()
Returns:

The matrix Q represented in terms of elementary reflectors.

inline HouseholderQ<T, L> thin_Q()
Returns:

The matrix Q1 ( \( Q = [Q_1 \: Q_2]\)) represented in terms of elementary reflectors.

inline TriangularMap<T, kUpper, L> R()
Returns:

the upper triangular matrix R.

inline TriangularMap<T, kUpper, L> thin_R()
Returns:

the upper triangular matrix R1 ( \( R = [R_1 \: 0]^T\)).

template<template<typename, typename> typename S>
inline void reorder(Matrix<T, L, S> &X)

Reorder the columns of the matrix X based on the permutation using during the column-pivoting QR factorization. Only available if kColPivot flag is set. Otherwise, this method result in a noop.

inline const mat_type &compressed_mat() const
Returns:

The matrix containing both the matrix R (on and above the diagonal) and the elementary reflectors that defines the matrix Q (below the diagonal).

inline const Vector<T, L> &reflect_coeff() const
Returns:

The coefficients of the elementary reflectors that defines the matrix Q.

inline const Vector<index_t, L> &pivots() const
Returns:

The index of the pivots using during the factorization. Only available if kColPivot flag is set.

inline index_t rows() const
Returns:

The number of rows of the underlying matrix

inline index_t cols() const
Returns:

The number of cols of the underlying matrix

LQ

template<Number T, DataLocation Loc, int Flags>
class LQ

A solver representing the LQ factorization of an \( m \times n \) matrix A:

\[ \mathbf{A} = \mathbf{LQ} \]
where Q is an orthogonal matrix and R is a lower triangular.

If \( m \leq n \), then the LQ factorization is given by

\[ \mathbf{A} = \begin{bmatrix} \mathbf{L} & 0 \end{bmatrix} \begin{bmatrix} \mathbf{Q_1} // \mathbf{Q_2} \end{bmatrix} \]
where L is an \( m \times m \) lower triangular matrix with real diagonal elements, and Q is an \( n \times n \) orthogonal (or unitary) matrix. The matrix \( \mathbf{L_1} \mathbf{Q_1} \) is referred as the thin LQ factorization of the matrix A.

Usage:

Matrix<double> A(m, n);
Matrix<double> B(m, n);
Matrix<double> Q(n, n);
LQ<double> F = lq(A);
Q = F.Q();
B = A * F.Q();

Note

Multiplying the matrix qr.thin_Q() directly as shown by the example is only allowed when m <= n. Otherwise, you need to form the matrix Q1 explicitly: Q1 = qr.thin_Q().

Note

When using column-pivoting (Flags & colPivot), the solution vector will not match the original problem. In this case, use the reorder() method to return the entries to their original positions.

Template Parameters:
  • T – Floating-point type

  • Loc – Location of the data.

  • Flags – Flags for the QR decomposition. Available: kDefaultFlags (default flags) and kInPlace (the factorization is done in-place, overwriting A).

Public Functions

inline LQ()

Default constructor.

template<template<typename, typename> typename S>
inline LQ(const Matrix<T, Loc, S> &A)

Creates a new solver and factorize the matrix A.

Parameters:

A[in] input matrix

virtual ~LQ() = default

Destructor.

template<template<typename, typename> typename S>
inline LQ &factorize(const Matrix<T, Loc, S> &A)

Calculates the LQ factorization of the matrix A.

Parameters:

A[in] input matrix

inline HouseholderQ<T, Loc> Q()
Returns:

The matrix Q represented in terms of elementary reflectors.

inline HouseholderQ<T, Loc> thin_Q()
Returns:

The matrix Q1 ( \( Q = [Q_1 \: Q_2]^T\)) represented in terms of elementary reflectors.

inline TriangularMap<T, kLower, Loc> L()
Returns:

the lower triangular matrix L.

inline TriangularMap<T, kLower, Loc> thin_L()
Returns:

the lower triangular matrix L1 ( \( L = [L_1 \: L_2]\)).

inline const mat_type &compressed_mat() const
Returns:

The matrix containing both the matrix R (on and above the diagonal) and the elementary reflectors that defines the matrix Q (below the diagonal).

inline const Vector<T, Loc> &reflect_coeff() const
Returns:

The coefficients of the elementary reflectors that defines the matrix Q.

inline index_t rows() const
Returns:

The number of rows of the underlying matrix

inline index_t cols() const
Returns:

The number of cols of the underlying matrix

Householder Q

template<Number ElemType, DataLocation Loc>
class HouseholderQ : public HighFM::internal::BLASOp<HouseholderQ<ElemType, Loc>>

The matrix Q from a QR or LQ decomposition represented in terms of Householder reflections.

Usage:

Matrix<float> A(m, n);
Matrix<float> Q(m, m);
Matrix<float> B(m, n);
Matrix<float> C(m, n);

auto F = qr(A);
Q = F.Q();              // Assemble the matrix Q
C = F.Q() * B;          // Apply the matrix Q without explicitly constructing it
C = F.Q().T() * B;      // Apply the transpose of the matrix Q
B = F.Q() * B;          // Apply the matrix Q in-place

Note

When using the thin (or partial) QR factorization, you must explicitly form the matrix Q before applying it to a matrix if Q is not square.

Template Parameters:
  • ElemType – Floating-point type

  • Loc – Location of the data.

Public Functions

template<int Flags>
inline HouseholderQ(const QR<ElemType, Loc, Flags> &fact, bool is_thin)

Creates a new HouseholderQ matrix from the QR factorization.

Parameters:
  • fact[in] An object containing the QR factorization

  • is_thin[in] Indicate if the object contains the full or the partial (or thin) QR factorization.

template<int Flags>
inline HouseholderQ(const LQ<ElemType, Loc, Flags> &fact, bool is_thin)

Creates a new HouseholderQ matrix from the QR factorization.

Parameters:
  • fact[in] An object containing the QR factorization

  • is_thin[in] Indicate if the object contains the full or the partial (or thin) QR factorization.

virtual ~HouseholderQ() = default

Destructor.

inline index_t rows() const noexcept
Returns:

The number of rows in the container.

inline index_t cols() const noexcept
Returns:

The number of columns in the container.

inline internal::ConjugateTransposeOp<HouseholderQ<ElemType, Loc>, internal::kTranspose> T()
Returns:

An abstract object representing the transpose of the matrix Q.

inline internal::ConjugateTransposeOp<HouseholderQ<ElemType, Loc>, internal::kAdjoint> H()
Returns:

An abstract object representing the adjoint of the matrix Q.