alpaqa 1.0.0a14
Nonconvex constrained optimization
Loading...
Searching...
No Matches
limited-memory-qr.hpp
Go to the documentation of this file.
1#pragma once
2
5#include <Eigen/Jacobi>
6#include <type_traits>
7
8namespace alpaqa {
9
10/// Incremental QR factorization using modified Gram-Schmidt with
11/// reorthogonalization.
12///
13/// Computes A = QR while allowing efficient removal of the first
14/// column of A or adding new columns at the end of A.
15template <Config Conf = DefaultConfig>
17 public:
19 LimitedMemoryQR() = default;
20
21 /// @param n
22 /// The size of the vectors, the number of rows of A.
23 /// @param m
24 /// The maximum number of columns of A.
25 ///
26 /// The maximum dimensions of Q are n×m and the maximum dimensions of R are
27 /// m×m.
29
30 length_t n() const { return Q.rows(); }
31 length_t m() const { return Q.cols(); }
32 length_t size() const { return n(); }
33 length_t history() const { return m(); }
34
35 /// Add the given column to the right.
36 template <class VecV>
37 void add_column(const VecV &v) {
38 assert(num_columns() < m());
39
40 auto q = Q.col(q_idx);
41 auto r = R.col(r_idx_end);
42
43 // Modified Gram-Schmidt to make q orthogonal to Q
44 q = v;
45 for (index_t i = 0; i < q_idx; ++i) {
46 real_t s = Q.col(i).dot(q);
47 r(i) = s;
48 q -= s * Q.col(i);
49 }
50
51 // Compute the norms of orthogonalized q and original v
52 real_t norm_q = q.norm();
53 real_t norm_v = v.norm();
54
55 // If ‖q‖ is significantly smaller than ‖v‖, perform
56 // reorthogonalization
57 auto η = real_t(0.7);
58 while (norm_q < η * norm_v) {
60 for (index_t i = 0; i < q_idx; ++i) {
61 real_t s = Q.col(i).dot(q);
62 r(i) += s;
63 q -= s * Q.col(i);
64 }
65 norm_v = norm_q;
66 norm_q = q.norm();
67 }
68
69 // Normalize q such that new matrix (Q q) remains orthogonal (i.e. has
70 // orthonormal columns)
71 r(q_idx) = norm_q;
72 q /= norm_q;
73 // Keep track of the minimum/maximum diagonal element of R
74 min_eig = std::min(min_eig, norm_q);
75 max_eig = std::max(max_eig, norm_q);
76
77 // Increment indices, add a column to Q and R.
78 ++q_idx;
80 }
81
82 /// Remove the leftmost column.
84 assert(num_columns() > 0);
85
86 // After removing the first column of the upper triangular matrix R,
87 // it becomes upper Hessenberg. Givens rotations are used to make it
88 // triangular again.
89 Eigen::JacobiRotation<real_t> G;
90 index_t r = 0; // row index of R
91 index_t c = r_succ(r_idx_start); // column index of R in storage
92 // Loop over the diagonal elements of R:
93 while (r < q_idx - 1) {
94 // Compute the Givens rotation that makes the subdiagonal element
95 // of column c of R zero.
96 G.makeGivens(R(r, c), R(r + 1, c), &R(r, c));
97 // Apply it to the remaining columns of R.
98 // Not the current column, because the diagonal element was updated
99 // by the makeGivens function, and the subdiagonal element doesn't
100 // have to be set to zero explicitly, it's implicit.
101 // Also not the previous columns, because they are already
102 // implicitly zero below the diagonal and this rotation wouldn't
103 // have any effect there.
104 // TODO: can this be sped up by applying it in two blocks instead
105 // of column by column?
106 for (index_t cc = r_succ(c); cc != r_idx_end; cc = r_succ(cc))
107 R.col(cc).applyOnTheLeft(r, r + 1, G.adjoint());
108 // Apply the inverse of the Givens rotation to Q.
109 Q.block(0, 0, Q.rows(), q_idx).applyOnTheRight(r, r + 1, G);
110 // Keep track of the minimum/maximum diagonal element of R
111 min_eig = std::min(min_eig, R(r, c));
112 max_eig = std::max(max_eig, R(r, c));
113 // Advance indices to next diagonal element of R.
114 ++r;
115 c = r_succ(c);
116 }
117 // Remove rightmost column of Q, since it corresponds to the bottom row
118 // of R, which was set to zero by the Givens rotations
119 --q_idx;
120 // Remove the first column of R.
122 }
123
124 /// Solve the least squares problem Ax = b.
125 /// Do not divide by elements that are smaller in absolute value than @p tol.
126 template <class VecB, class VecX>
127 void solve_col(const VecB &b, VecX &x, real_t tol = 0) const {
128 // Iterate over the diagonal of R, starting at the bottom right,
129 // this is standard back substitution
130 // (recall that R is stored in a circular buffer, so R.col(i) is
131 // not the mathematical i-th column)
132 auto rev_bgn = ring_reverse_iter().begin();
133 auto rev_end = ring_reverse_iter().end();
134 auto fwd_end = ring_iter().end();
135 for (auto it_d = rev_bgn; it_d != rev_end; ++it_d) {
136 // Row/column index of diagonal element of R
137 auto [rR, cR] = *it_d;
138 // Don't divide by very small diagonal elements
139 if (std::abs(R(rR, cR)) < tol) {
140 x(rR) = real_t{0};
141 continue;
142 }
143 // (r is the zero-based mathematical index, c is the index in
144 // the circular buffer)
145 x(rR) = Q.col(rR).transpose() * b; // Compute rhs Qᵀb
146 // In the current row of R, iterate over the elements to the
147 // right of the diagonal
148 // Iterating from left to right seems to give better results
149 for (auto it_c = it_d.forwardit; it_c != fwd_end; ++it_c) {
150 auto [rX2, cR2] = *it_c;
151 x(rR) -= R(rR, cR2) * x(rX2);
152 }
153 x(rR) /= R(rR, cR); // Divide by diagonal element
154 }
155 }
156
157 /// Solve the least squares problem AX = B.
158 /// Do not divide by elements that are smaller in absolute value than @p tol.
159 template <class MatB, class MatX>
160 void solve(const MatB &B, MatX &X, real_t tol = 0) const {
161 assert(B.cols() <= X.cols());
162 assert(B.rows() >= Q.rows());
163 assert(X.rows() >= Eigen::Index(num_columns()));
164 // Each column of the right hand side is solved as an individual system
165 for (Eigen::Index cB = 0; cB < B.cols(); ++cB) {
166 auto b = B.col(cB);
167 auto x = X.col(cB);
168 solve_col(b, x, tol);
169 }
170 }
171
172 template <class Derived>
173 using solve_ret_t = std::conditional_t<
174 Eigen::internal::traits<Derived>::ColsAtCompileTime == 1, vec, mat>;
175
176 /// Solve the least squares problem AX = B.
177 template <class Derived>
178 solve_ret_t<Derived> solve(const Eigen::DenseBase<Derived> &B) {
179 solve_ret_t<Derived> X(m(), B.cols());
180 solve(B, X);
181 return X;
182 }
183
184 /// Get the full, raw storage for the orthogonal matrix Q.
185 const mat &get_raw_Q() const { return Q; }
186 /// Get the full, raw storage for the upper triangular matrix R.
187 /// The columns of this matrix are permuted because it's stored as a
188 /// circular buffer for efficiently appending columns to the end and
189 /// popping columns from the front.
190 const mat &get_raw_R() const { return R; }
191
192 /// Get the full storage for the upper triangular matrix R but with the
193 /// columns in the correct order.
194 /// @note Meant for tests only, creates a permuted copy.
195 mat get_full_R() const {
196 if (r_idx_start == 0)
197 return R;
198 // Using a permutation matrix here isn't as efficient as rotating the
199 // matrix manually, but this function is only used in tests, so it
200 // shouldn't matter.
201 Eigen::PermutationMatrix<Eigen::Dynamic> P(R.cols());
202 P.setIdentity();
203 std::rotate(P.indices().data(), P.indices().data() + r_idx_start,
204 P.indices().data() + P.size());
205 return R * P;
206 }
207 /// Get the matrix R such that Q times R is the original matrix.
208 /// @note Meant for tests only, creates a permuted copy.
209 mat get_R() const {
210 return get_full_R()
211 .block(0, 0, q_idx, q_idx)
213 }
214 /// Get the matrix Q such that Q times R is the original matrix.
215 /// @note Meant for tests only, creates a copy.
216 mat get_Q() const { return Q.block(0, 0, n(), q_idx); }
217
218 /// Multiply the matrix R by a scalar.
220 for (auto [i, r_idx] : ring_iter())
221 R.col(r_idx).topRows(i + 1) *= scal;
222 min_eig *= scal;
223 max_eig *= scal;
224 }
225
226 /// Get the number of MGS reorthogonalizations.
227 unsigned long get_reorth_count() const { return reorth_count; }
228 /// Reset the number of MGS reorthogonalizations.
230
231 /// Get the minimum eigenvalue of R.
232 real_t get_min_eig() const { return min_eig; }
233 /// Get the maximum eigenvalue of R.
234 real_t get_max_eig() const { return max_eig; }
235
236 /// Reset all indices, clearing the Q and R matrices.
237 void reset() {
238 q_idx = 0;
239 r_idx_start = 0;
240 r_idx_end = 0;
241 reorth_count = 0;
244 }
245
246 /// Re-allocate storage for a problem with a different size. Causes
247 /// a @ref reset.
249 Q.resize(n, m);
250 R.resize(m, m);
251 reset();
252 }
253
254 /// Get the number of columns that are currently stored.
255 length_t num_columns() const { return q_idx; }
256 /// Get the head index of the circular buffer (points to the oldest
257 /// element).
258 index_t ring_head() const { return r_idx_start; }
259 /// Get the tail index of the circular buffer (points to one past the most
260 /// recent element).
261 index_t ring_tail() const { return r_idx_end; }
262 /// Get the next index in the circular buffer.
263 index_t ring_next(index_t i) const { return r_succ(i); }
264 /// Get the previous index in the circular buffer.
265 index_t ring_prev(index_t i) const { return r_pred(i); }
266 /// Get the number of columns currently stored in the buffer.
267 length_t current_history() const { return q_idx; }
268
269 /// Get iterators in the circular buffer.
271 return {q_idx, r_idx_start, r_idx_end, m()};
272 }
273 /// Get reverse iterators in the circular buffer.
277
278 private:
279 mat Q; ///< Storage for orthogonal factor Q.
280 mat R; ///< Storage for upper triangular factor R.
281
282 index_t q_idx = 0; ///< Number of columns of Q being stored.
283 index_t r_idx_start = 0; ///< Index of the first column of R.
284 index_t r_idx_end = 0; ///< Index of the one-past-last column of R.
285
286 unsigned long reorth_count = 0; ///< Number of MGS reorthogonalizations.
287
288 real_t min_eig = +inf<config_t>; ///< Minimum eigenvalue of R.
289 real_t max_eig = -inf<config_t>; ///< Maximum eigenvalue of R.
290
291 /// Get the next index in the circular storage for R.
292 index_t r_succ(index_t i) const { return i + 1 < m() ? i + 1 : 0; }
293 /// Get the previous index in the circular storage for R.
294 index_t r_pred(index_t i) const { return i == 0 ? m() - 1 : i - 1; }
295};
296
297} // namespace alpaqa
Incremental QR factorization using modified Gram-Schmidt with reorthogonalization.
mat R
Storage for upper triangular factor R.
index_t r_idx_start
Index of the first column of R.
real_t max_eig
Maximum eigenvalue of R.
CircularRange< index_t > ring_iter() const
Get iterators in the circular buffer.
mat Q
Storage for orthogonal factor Q.
LimitedMemoryQR(length_t n, length_t m)
length_t num_columns() const
Get the number of columns that are currently stored.
mat get_Q() const
Get the matrix Q such that Q times R is the original matrix.
void solve(const MatB &B, MatX &X, real_t tol=0) const
Solve the least squares problem AX = B.
length_t current_history() const
Get the number of columns currently stored in the buffer.
unsigned long reorth_count
Number of MGS reorthogonalizations.
void remove_column()
Remove the leftmost column.
mat get_R() const
Get the matrix R such that Q times R is the original matrix.
index_t q_idx
Number of columns of Q being stored.
index_t r_pred(index_t i) const
Get the previous index in the circular storage for R.
void add_column(const VecV &v)
Add the given column to the right.
void solve_col(const VecB &b, VecX &x, real_t tol=0) const
Solve the least squares problem Ax = b.
solve_ret_t< Derived > solve(const Eigen::DenseBase< Derived > &B)
Solve the least squares problem AX = B.
real_t min_eig
Minimum eigenvalue of R.
void clear_reorth_count()
Reset the number of MGS reorthogonalizations.
index_t ring_tail() const
Get the tail index of the circular buffer (points to one past the most recent element).
real_t get_max_eig() const
Get the maximum eigenvalue of R.
index_t r_succ(index_t i) const
Get the next index in the circular storage for R.
std::conditional_t< Eigen::internal::traits< Derived >::ColsAtCompileTime==1, vec, mat > solve_ret_t
index_t ring_prev(index_t i) const
Get the previous index in the circular buffer.
index_t ring_head() const
Get the head index of the circular buffer (points to the oldest element).
const mat & get_raw_Q() const
Get the full, raw storage for the orthogonal matrix Q.
const mat & get_raw_R() const
Get the full, raw storage for the upper triangular matrix R.
void resize(length_t n, length_t m)
Re-allocate storage for a problem with a different size.
void reset()
Reset all indices, clearing the Q and R matrices.
index_t ring_next(index_t i) const
Get the next index in the circular buffer.
unsigned long get_reorth_count() const
Get the number of MGS reorthogonalizations.
real_t get_min_eig() const
Get the minimum eigenvalue of R.
void scale_R(real_t scal)
Multiply the matrix R by a scalar.
ReverseCircularRange< index_t > ring_reverse_iter() const
Get reverse iterators in the circular buffer.
mat get_full_R() const
Get the full storage for the upper triangular matrix R but with the columns in the correct order.
index_t r_idx_end
Index of the one-past-last column of R.
#define USING_ALPAQA_CONFIG(Conf)
Definition config.hpp:56
typename Conf::mat mat
Definition config.hpp:71
typename Conf::real_t real_t
Definition config.hpp:65
typename Conf::index_t index_t
Definition config.hpp:77
typename Conf::length_t length_t
Definition config.hpp:76
constexpr const auto inf
Definition config.hpp:85
typename Conf::vec vec
Definition config.hpp:66