Linear Algebra  arduino
Accessible implementations of linear algebra algorithms
HouseholderQR.cpp
Go to the documentation of this file.
1 #ifndef ARDUINO
3 #else
5 #endif
6 
7 #include <cassert>
8 
27  // For the intermediate calculations, we'll be working with RW.
28  // It is initialized to the rectangular matrix to be factored.
29  // At the end of this function, RW will contain the strict
30  // upper-triangular part of the matrix R (without the diagonal),
31  // and the complete scaled matrix of reflection vectors W, which is a
32  // lower-triangular matrix. The diagonal of R is stored separately in
33  // R_diag.
34 
35  assert(RW.rows() >= RW.cols());
36  assert(R_diag.size() == RW.cols());
37 
38  // Helper function to square a number
39  auto sq = [](double x) { return x * x; };
40 
41  for (size_t k = 0; k < RW.cols(); ++k) {
42  // Introduce a column vector x = A[k:M,k], it's the lower part of the
43  // k-th column of the matrix.
44  // First compute the norm of x:
45 
46  double sq_norm_x = 0;
47  for (size_t i = k; i < RW.rows(); ++i)
48  sq_norm_x += sq(RW(i, k));
49  double norm_x = std::sqrt(sq_norm_x);
50 
51  // x consists of two parts: its first element, x₀, and the rest, xₛ
52  // x = (x₀, xₛ)
53  // You can express the norm of x in terms of the norms of the two parts:
54  // ‖x‖² = x₀² + ‖xₛ‖²
55  double &x_0 = RW(k, k);
56 
57  // The goal of QR factorization is to introduce zeros below the diagonal
58  // in the R factor by transforming the vector x to a new vector that is
59  // all zero, except for the first component.
60  // Call this vector xₕ, the Householder reflection of x.
61  // Since the transformation has to be unitary (Q is a unitary matrix),
62  // this operation has to preserve the 2-norm. This means that the
63  // nonzero component of xₕ has to have the same 2-norm (energy) as x:
64  // xₕ = (±‖x‖, 0, ..., 0) = ±‖x‖·̅e₁
65  // where ̅e₁ is the first standard basis vector (1, 0, ..., 0).
66  //
67  // There are two (real) vectors that have the same energy as x in their
68  // first component, because the sign doesn't affect the energy.
69  // For numerical reasons, it's best to pick the sign opposite to the
70  // sign of the first component of x, so that x and xₕ are far apart.
71  //
72  // The reflector vector vₖ is the difference between the original
73  // vector x and its Householder reflection xₕ:
74  // vₖ = x - xₕ
75  // = x - (-sign(x₀)·‖x‖·̅e₁)
76  // = x + sign(x₀)·‖x‖·̅e₁
77  // = (x₀ + sign(x₀)·‖x‖, xₛ)
78  //
79  // Since vₖ will later be used to construct a projection matrix, it
80  // should be normalized.
81  // Computing the full norm is not necessary, since
82  // ‖vₖ‖² = (x₀ + sign(x₀)·‖x‖)² + ‖xₛ‖²
83  // = x₀² + 2·x₀·sign(x₀)·‖x‖ + ‖x‖² + ‖xₛ‖²
84  // = 2·x₀·sign(x₀)·‖x‖ + 2·‖x‖²
85  // = 2·|x₀|·‖x‖ + 2·‖x‖²
86  // ‖vₖ‖ = √2·√(|x₀|·‖x‖ + ‖x‖²)
87  //
88  // Normalize vₖ and call this uₖ:
89  // uₖ = vₖ / ‖vₖ‖
90  //
91  // For reasons that will become apparent in a moment, we'll keep the
92  // factor of √2 in our normalization of vₖ, call this vector wₖ:
93  // wₖ = √2·uₖ
94  // = √2·vₖ / ‖vₖ‖
95  // = vₖ / √(|x₀|·‖x‖ + ‖x‖²)
96  // Note how the sum only adds up numbers with the same sign. This
97  // prevents catastrophic cancelations.
98  //
99  // It is clear that if ‖x‖ = 0, this normalization will fail. In that
100  // case, we set uₖ = ̅e₁ or wₖ = √2·̅e₁.
101  //
102  // x will be overwritten by wₖ. The vector xₕ only has a single nonzero
103  // component. It is saved in the R_diag vector.
104 
105  if (norm_x >= std::numeric_limits<double>::min() * 2) {
106  double x_p = -std::copysign(norm_x, x_0); // -sign(x₀)·‖x‖
107  double v_0 = x_0 - x_p;
108  double norm_v_sq2 = std::sqrt(std::abs(x_0) * norm_x + sq_norm_x);
109 
110  // Overwrite x with vₖ:
111  x_0 = v_0;
112  // the other components of x (xₛ) are already equal to the bottom
113  // part of vₖ, so they don't have to be overwritten explicitly.
114 
115  // Then normalize x (= vₖ) to obtain wₖ:
116  for (size_t i = k; i < RW.rows(); ++i)
117  RW(i, k) /= norm_v_sq2;
118 
119  // Save the first component of xₕ:
120  R_diag(k) = x_p;
121  } else {
122  // Overwrite x with wₖ = √2·̅e₁:
123  x_0 = std::sqrt(2);
124  // the other components of x (xₛ) are already equal to zero, since
125  // ‖x‖ = 0.
126 
127  // Save the first component of xₕ:
128  R_diag(k) = 0;
129  }
130 
131  // Now that the reflection vector vₖ (wₖ) is known, the rest of the
132  // matrix A can be updated, that is A[k:m,k:n].
133  //
134  // The reflection was defined in terms of the reflection vector
135  // vₖ = x - xₕ. To reflect x onto xₕ, you can first project x onto the
136  // orthogonal complement of the space spanned by vₖ, call this
137  // projection xₚ, and then add the difference between xₚ and x to xₚ
138  // again, in order to reflect about the orthogonal complement of vₖ,
139  // rather than projecting onto it:
140  // xₕ = xₚ + (xₚ - x)
141  // = 2·xₚ - x
142  // The projection matrix onto the orthogonal complement of the space
143  // spanned by vₖ is:
144  // P⟂ = I - uₖ·uₖᵀ
145  // where I is the identity matrix, and uₖ the unit vector in the
146  // direction of vₖ, as defined above.
147  // This allows us to compute xₚ:
148  // xₚ = P⟂·x
149  // = I·x - uₖ·uₖᵀ·x
150  // = x - uₖ·uₖᵀ·x
151  // xₕ = 2·xₚ - x
152  // = 2·x - 2·uₖ·uₖᵀ·x - x
153  // = x - 2·uₖ·uₖᵀ·x
154  // Because of our choice of wₖ = √2·uₖ earlier, this simplifies to
155  // xₕ = x - wₖ·wₖᵀ·x
156  // The Householder reflector is thus the linear transformation given by
157  // H = I - wₖ·wₖᵀ
158  //
159  // The final step is to apply this reflector to the remaining part of
160  // the matrix A, i.e. A[k:m,k:n].
161  // The first column x = A[k:m,k] has already been updated, that's how
162  // the reflector was computed, so only the submatrix A[k+1:m,k:n] has
163  // to be updated.
164  // A'[k+1:m,k:n] = H·A[k+1:m,k:n]
165  // = A[k+1:m,k:n] - wₖ·wₖᵀ·A[k+1:m,k:n]
166  // This can be computed column-wise:
167  // aᵢ' = aᵢ - wₖ·wₖᵀ·aᵢ
168  // where aᵢ is the i-th column of A.
169 
170  for (size_t c = k + 1; c < RW.cols(); ++c) {
171  // Compute wₖᵀ·aᵢ
172  double dot_product = 0;
173  for (size_t r = k; r < RW.rows(); ++r)
174  dot_product += RW(r, k) * RW(r, c);
175  // Subtract wₖ·wₖᵀ·aᵢ
176  for (size_t r = k; r < RW.rows(); ++r)
177  RW(r, c) -= RW(r, k) * dot_product;
178  }
179  }
180  state = Factored;
181 }
183 
190  assert(is_factored());
191  assert(RW.rows() == B.rows());
192  // Apply the Householder reflectors to each column of B.
193  for (size_t i = 0; i < B.cols(); ++i) {
194  // Recall that the Householder reflector H is applied as follows:
195  // bᵢ'[k+1:m] = H·bᵢ[k+1:m]
196  // = bᵢ[k+1:m] - wₖ·wₖᵀ·bᵢ[k+1:m]
197  for (size_t k = 0; k < RW.cols(); ++k) {
198  // Compute wₖᵀ·bᵢ
199  double dot_product = 0;
200  for (size_t r = k; r < RW.rows(); ++r)
201  dot_product += RW(r, k) * B(r, i);
202  // Subtract wₖ·wₖᵀ·bᵢ
203  for (size_t r = k; r < RW.rows(); ++r)
204  B(r, i) -= RW(r, k) * dot_product;
205  }
206  }
207 }
209 
216  assert(is_factored());
217  assert(RW.rows() == X.rows());
218  // Apply the Householder reflectors in reverse order to each column of X.
219  for (size_t i = 0; i < X.cols(); ++i) {
220  // Recall that the Householder reflector H is applied as follows:
221  // xᵢ'[k+1:m] = H·xᵢ[k+1:m]
222  // = xᵢ[k+1:m] - wₖ·wₖᵀ·xᵢ[k+1:m]
223  for (size_t k = RW.cols(); k-- > 0;) {
224  // Compute wₖᵀ·xᵢ
225  double dot_product = 0;
226  for (size_t r = k; r < RW.rows(); ++r)
227  dot_product += RW(r, k) * X(r, i);
228  // Subtract wₖ·wₖᵀ·xᵢ
229  for (size_t r = k; r < RW.rows(); ++r)
230  X(r, i) -= RW(r, k) * dot_product;
231  }
232  }
233 }
235 
241 void HouseholderQR::back_subs(const Matrix &B, Matrix &X) const {
242  // Solve upper triangular system RX = B by solving each column of B as a
243  // vector system Rxᵢ = bᵢ
244  //
245  // ┌ ┐┌ ┐ ┌ ┐
246  // │ r₁₁ r₁₂ r₁₃ r₁₄ ││ x₁ᵢ │ │ b₁ᵢ │
247  // │ r₂₂ r₂₃ r₂₄ ││ x₂ᵢ │ = │ b₂ᵢ │
248  // │ r₃₃ r₃₄ ││ x₃ᵢ │ │ b₃ᵢ │
249  // │ r₄₄ ││ x₄ᵢ │ │ b₄ᵢ │
250  // └ ┘└ ┘ └ ┘
251  //
252  // b₄ᵢ = r₄₄·x₄ᵢ ⟺ x₄ᵢ = b₄ᵢ/r₄₄
253  // b₃ᵢ = r₃₃·x₃ᵢ + r₃₄·x₄ᵢ ⟺ x₃ᵢ = (b₃ᵢ - r₃₄·x₄ᵢ)/r₃₃
254  // b₂ᵢ = r₂₂·x₂ᵢ + r₂₃·x₃ᵢ + r₂₄·x₄ᵢ ⟺ x₂ᵢ = (b₂ᵢ - r₂₃·x₃ᵢ + r₂₄·x₄ᵢ)/r₂₂
255  // ...
256 
257  for (size_t i = 0; i < B.cols(); ++i) {
258  for (size_t k = RW.cols(); k-- > 0;) {
259  X(k, i) = B(k, i);
260  for (size_t j = k + 1; j < RW.cols(); ++j) {
261  X(k, i) -= RW(k, j) * X(j, i);
262  }
263  X(k, i) /= R_diag(k);
264  }
265  }
266 }
268 
275  // If AX = B, then QRX = B, or RX = QᵀB, so first apply Qᵀ to B:
276  apply_QT_inplace(B);
277 
278  // To solve RX = QᵀB, use back substitution:
279 
280  // If the matrix is square, B and X have the same size, so we can reuse B's
281  // storage if we operate on B directly:
282  if (RW.cols() == RW.rows()) {
283  back_subs(B, B);
284  }
285  // If the matrix is rectangular, the sizes of B and X differ, so use a
286  // separate result variable:
287  else {
288  Matrix X(RW.cols(), B.cols());
289  back_subs(B, X);
290  B = std::move(X);
291  }
292 }
294 
295 // All implementations of the less interesting functions can be found here:
Vector R_diag
Contains the diagonal elements of R.
void apply_QT_inplace(Matrix &B) const
Compute the product QᵀB, overwriting B with the result.
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...
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.
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
size_t cols() const
Get the number of columns of the matrix.
Definition: Matrix.hpp:83
size_t size() const
Get the number of elements in the vector.
Definition: Matrix.hpp:310