Linear Algebra  arduino
Accessible implementations of linear algebra algorithms
HouseholderQR.ipp
Go to the documentation of this file.
1 #ifndef ARDUINO
3 #else
5 #endif
6 
7 #include <cassert>
8 
9 void HouseholderQR::compute(Matrix &&matrix) {
10  RW = std::move(matrix);
11  R_diag.resize(RW.cols());
13 }
14 
15 void HouseholderQR::compute(const Matrix &matrix) {
16  RW = matrix;
17  R_diag.resize(RW.cols());
19 }
20 
22  Matrix result = B;
23  apply_QT_inplace(result);
24  return result;
25 }
26 
29  return std::move(B);
30 }
31 
33  Matrix result = X;
34  apply_Q_inplace(result);
35  return result;
36 }
37 
39  apply_Q_inplace(B);
40  return std::move(B);
41 }
42 
44  assert(is_factored());
45  assert(R.rows() == RW.rows());
46  assert(R.cols() == RW.cols());
47  for (size_t r = 0; r < R.cols(); ++r) {
48  for (size_t c = 0; c < r; ++c)
49  R(r, c) = 0;
50  R(r, r) = R_diag(r);
51  for (size_t c = r + 1; c < R.cols(); ++c)
52  R(r, c) = RW(r, c);
53  }
54  for (size_t r = R.cols(); r < R.rows(); ++r) {
55  for (size_t c = 0; c < R.cols(); ++c)
56  R(r, c) = 0;
57  }
58 }
59 
61  Matrix R(RW.rows(), RW.cols());
62  get_R_inplace(R);
63  return R;
64 }
65 
68  for (size_t r = 0; r < RW.cols(); ++r) {
69  for (size_t c = 0; c < r; ++c)
70  RW(r, c) = 0;
71  RW(r, r) = R_diag(r);
72  }
73  for (size_t r = RW.cols(); r < RW.rows(); ++r) {
74  for (size_t c = 0; c < RW.cols(); ++c)
75  RW(r, c) = 0;
76  }
77  return std::move(RW);
78 }
79 
81  assert(Q.rows() == RW.rows());
82  assert(Q.cols() == RW.rows());
83  Q.fill_identity();
84  apply_Q_inplace(Q);
85 }
86 
88  SquareMatrix Q(RW.rows());
89  get_Q_inplace(Q);
90  return Q;
91 }
92 
94  Matrix B_cpy = apply_QT(B);
95  Matrix X(RW.cols(), B.cols());
96  back_subs(B_cpy, X);
97  return X;
98 }
99 
101  solve_inplace(B);
102  return std::move(B);
103 }
104 
106  return Vector(solve(static_cast<const Matrix &>(b)));
107 }
108 
110  solve_inplace(b);
111  return std::move(b);
112 }
113 
114 #ifndef NO_IOSTREAM_SUPPORT
115 
116 #include <iomanip>
117 #include <iostream>
118 
119 // LCOV_EXCL_START
120 
121 std::ostream &operator<<(std::ostream &os, const HouseholderQR &qr) {
122  if (!qr.is_factored()) {
123  os << "Not factored." << std::endl;
124  return os;
125  }
126 
127  Matrix Q = qr.get_Q();
128  os << "Q = " << std::endl;
129  Q.print(os);
130 
131  // Output field width (characters)
132  int w = os.precision() + 9;
133 
134  const auto &RW = qr.get_RW();
135  const auto &R_diag = qr.get_R_diag();
136  os << "R = " << std::endl;
137  for (size_t r = 0; r < RW.cols(); ++r) {
138  for (size_t c = 0; c < r; ++c)
139  os << std::setw(w) << 0;
140  os << std::setw(w) << R_diag(r);
141  for (size_t c = r + 1; c < RW.cols(); ++c)
142  os << std::setw(w) << RW(r, c);
143  os << std::endl;
144  }
145  for (size_t r = RW.cols(); r < RW.rows(); ++r) {
146  for (size_t c = 0; c < RW.cols(); ++c)
147  os << std::setw(w) << 0;
148  }
149  return os;
150 }
151 
152 // LCOV_EXCL_STOP
153 
154 #endif
QR factorization using Householder reflectors.
Vector R_diag
Contains the diagonal elements of R.
std::ostream & operator<<(std::ostream &os, const HouseholderQR &qr)
Print the Q and R matrices of a HouseholderQR object.
void apply_QT_inplace(Matrix &B) const
Compute the product QᵀB, overwriting B with the result.
void get_Q_inplace(SquareMatrix &Q) const
Compute the unitary matrix Q and copy it to the given matrix.
Matrix solve(const Matrix &B) const
Solve the system AX = B or QRX = B.
Matrix apply_Q(const Matrix &X) const
Compute the product QB.
Matrix get_R() const &
Get a copy of the upper-triangular matrix R.
const Vector & get_R_diag() const &
Get the internal storage of the diagonal elements of R.
Matrix && steal_R()
Get the upper-triangular matrix R, reusing the internal storage.
enum HouseholderQR::State state
void back_subs(const Matrix &B, Matrix &X) const
Back substitution algorithm for solving upper-triangular systems RX = B.
Matrix RW
Result of a Householder QR factorization: stores the strict upper-triangular part of matrix R and the...
Matrix apply_QT(const Matrix &B) const
Compute the product QᵀB.
void solve_inplace(Matrix &B) const
Solve the system AX = B or QRX = B.
bool is_factored() const
Check if this object contains a valid factorization.
void compute_factorization()
The actual QR factorization algorithm.
SquareMatrix get_Q() const
Compute the unitary matrix Q.
const Matrix & get_RW() const &
Get the internal storage of the strict upper-triangular part of R and the Householder reflector vecto...
void compute(Matrix &&matrix)
Perform the QR factorization of the given matrix.
void get_R_inplace(Matrix &R) const
Copy the upper-triangular matrix R to the given matrix.
void apply_Q_inplace(Matrix &X) const
Compute the product QB, overwriting B with the result.
General matrix class.
Definition: Matrix.hpp:34
size_t rows() const
Get the number of rows of the matrix.
Definition: Matrix.hpp:81
void fill_identity()
Fill the matrix as an identity matrix (all zeros except the diagonal which is one).
Definition: Matrix.cpp:129
size_t cols() const
Get the number of columns of the matrix.
Definition: Matrix.hpp:83
void print(std::ostream &os, uint8_t precision=0, uint8_t width=0) const
Print a matrix.
Definition: Matrix.cpp:242
Square matrix class.
Definition: Matrix.hpp:572
A column vector (n×1 matrix).
Definition: Matrix.hpp:278
void resize(size_t size)
Resize the vector.
Definition: Matrix.hpp:307