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