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