PANOC-ALM  quadratic-penalty
Nonconvex constrained optimization
limited-memory-qr.hpp
Go to the documentation of this file.
1 #include <Eigen/Jacobi>
2 #include <cstddef>
4 #include <panoc-alm/util/vec.hpp>
5 #include <type_traits>
6 
7 namespace pa {
8 
15  public:
16  LimitedMemoryQR() = default;
17 
25  LimitedMemoryQR(size_t n, size_t m) : Q(n, m), R(m, m) {}
26 
27  size_t n() const { return Q.rows(); }
28  size_t m() const { return Q.cols(); }
29 
31  template <class VecV>
32  void add_column(const VecV &v) {
33  assert(q_idx < m());
34 
35  auto q = Q.col(q_idx);
36  auto r = R.col(r_idx_end);
37 
38  // Modified Gram-Schmidt to make q orthogonal to Q
39  q = v;
40  for (size_t i = 0; i < q_idx; ++i) {
41  real_t s = Q.col(i).dot(q);
42  r(i) = s;
43  q -= s * Q.col(i);
44  }
45 
46  // Compute the norms of orthogonalized q and original v
47  real_t norm_q = q.norm();
48  real_t norm_v = v.norm();
49 
50  // If ‖q‖ is significantly smaller than ‖v‖, perform
51  // reorthogonalization
52  real_t η = 0.7;
53  while (norm_q < η * norm_v) {
54  ++reorth_count;
55  for (size_t i = 0; i < q_idx; ++i) {
56  real_t s = Q.col(i).dot(q);
57  r(i) += s;
58  q -= s * Q.col(i);
59  }
60  norm_v = norm_q;
61  norm_q = q.norm();
62  }
63 
64  // Normalize q such that new matrix (Q q) remains orthogonal (i.e. has
65  // orthonormal columns)
66  r(q_idx) = norm_q;
67  q /= norm_q;
68 
69  // Increment indices, add a column to Q and R.
70  ++q_idx;
72  }
73 
75  void remove_column() {
76  assert(q_idx > 0);
77 
78  // After removing the first colomn of the upper triangular matrix R,
79  // it becomes upper Hessenberg. Givens rotations are used to make it
80  // triangular again.
81  Eigen::JacobiRotation<real_t> G;
82  size_t r = 0; // row index of R
83  size_t c = r_succ(r_idx_start); // column index of R in storage
84  while (r < q_idx - 1) {
85  // Compute the Givens rotation that makes the subdiagonal element
86  // of column c or R zero.
87  G.makeGivens(R(r, c), R(r + 1, c), &R(r, c));
88  // Apply it to the remaining columns of R.
89  // Not the current column, because the diagonal element was updated
90  // by the makeGivens function, and the subdiagonal element doesn't
91  // have to be set to zero explicitly, it's implicit.
92  // Also not the previous columns, because they are already
93  // implicitly zero below the diagonal and this rotation wouldn't
94  // have any effect there.
95  // TODO: can this be sped up by applying it in two blocks instead
96  // of column by column?
97  for (size_t cc = r_succ(c); cc != r_idx_end; cc = r_succ(cc))
98  R.col(cc).applyOnTheLeft(r, r + 1, G.adjoint());
99  // Apply the inverse of the Givens rotation to Q.
100  Q.block(0, 0, Q.rows(), q_idx).applyOnTheRight(r, r + 1, G);
101  // Advance indices to next diagonal element of R.
102  ++r;
103  c = r_succ(c);
104  }
105  // Remove rightmost column of Q, since it corresponds to the bottom row
106  // of R, which was set to zero by the Givens rotations
107  --q_idx;
108  // Remove the first column of R.
110  }
111 
113  template <class VecB, class VecX>
114  void solve_col(const VecB &b, VecX &x) const {
115  // Iterate over the diagonal of R, starting at the bottom right,
116  // this is standard back substitution
117  // (recall that R is stored in a circular buffer, so R.col(i) is
118  // not the mathematical i-th column)
119  auto rev_bgn = ring_reverse_iter().begin();
120  auto rev_end = ring_reverse_iter().end();
121  auto fwd_end = ring_iter().end();
122  for (auto it_d = rev_bgn; it_d != rev_end; ++it_d) {
123  // Row/column index of diagonal element of R
124  auto [rR, cR] = *it_d;
125  // (r is the zero-based mathematical index, c is the index in
126  // the circular buffer)
127  x(rR) = Q.col(rR).transpose() * b; // Compute rhs Qᵀb
128  // In the current row of R, iterate over the elements to the
129  // right of the diagonal
130  // Iterating from left to right seems to give better results
131  for (auto it_c = it_d.forwardit; it_c != fwd_end; ++it_c) {
132  auto [rX2, cR2] = *it_c;
133  x(rR) -= R(rR, cR2) * x(rX2);
134  }
135  x(rR) /= R(rR, cR); // Divide by diagonal element
136  }
137  }
138 
140  template <class MatB, class MatX>
141  void solve(const MatB &B, MatX &X) const {
142  assert(B.cols() <= X.cols());
143  assert(B.rows() >= Q.rows());
144  assert(X.rows() >= Eigen::Index(num_columns()));
145  // Each column of the right hand side is solved as an individual system
146  for (Eigen::Index cB = 0; cB < B.cols(); ++cB) {
147  auto b = B.col(cB);
148  auto x = X.col(cB);
149  solve_col(b, x);
150  }
151  }
152 
153  template <class Derived>
154  using solve_ret_t = std::conditional_t<
155  Eigen::internal::traits<Derived>::ColsAtCompileTime == 1, vec, mat>;
156 
158  template <class Derived>
159  solve_ret_t<Derived> solve(const Eigen::DenseBase<Derived> &B) {
160  solve_ret_t<Derived> X(m(), B.cols());
161  solve(B, X);
162  return X;
163  }
164 
166  const mat &get_raw_Q() const { return Q; }
171  const mat &get_raw_R() const { return R; }
172 
176  mat get_full_R() const {
177  if (r_idx_start == 0)
178  return R;
179  // Using a permutation matrix here isn't as efficient as rotating the
180  // matrix manually, but this function is only used in tests, so it
181  // shouldn't matter.
182  Eigen::PermutationMatrix<Eigen::Dynamic> P(R.cols());
183  P.setIdentity();
184  std::rotate(P.indices().data(), P.indices().data() + r_idx_start,
185  P.indices().data() + P.size());
186  return R * P;
187  }
190  mat get_R() const {
191  return get_full_R()
192  .block(0, 0, q_idx, q_idx)
193  .triangularView<Eigen::Upper>();
194  }
197  mat get_Q() const { return Q.block(0, 0, n(), q_idx); }
198 
200  void scale_R(real_t scal) {
201  for (auto [i, r_idx] : ring_iter())
202  R.col(r_idx).topRows(i + 1) *= scal;
203  }
204 
206  size_t get_reorth_count() const { return reorth_count; }
209 
211  void reset() {
212  q_idx = 0;
213  r_idx_start = 0;
214  r_idx_end = 0;
215  reorth_count = 0;
216  }
217 
220  void resize(size_t n, size_t m) {
221  Q.resize(n, m);
222  R.resize(m, m);
223  reset();
224  }
225 
227  size_t num_columns() const { return q_idx; }
230  size_t ring_head() const { return r_idx_start; }
233  size_t ring_tail() const { return r_idx_end; }
235  size_t ring_next(size_t i) const { return r_succ(i); }
237  size_t ring_prev(size_t i) const { return r_pred(i); }
238 
241  return {q_idx, r_idx_start, r_idx_end, m()};
242  }
245  return ring_iter();
246  }
247 
248  private:
249  mat Q;
250  mat R;
251 
252  size_t q_idx = 0;
253  size_t r_idx_start = 0;
254  size_t r_idx_end = 0;
255 
256  size_t reorth_count = 0;
257 
259  size_t r_succ(size_t i) const { return i + 1 < m() ? i + 1 : 0; }
261  size_t r_pred(size_t i) const { return i == 0 ? m() - 1 : i - 1; }
262 };
263 
264 } // namespace pa
pa::LimitedMemoryQR::m
size_t m() const
Definition: limited-memory-qr.hpp:28
pa::LimitedMemoryQR::Q
mat Q
Storage for orthogonal factor Q.
Definition: limited-memory-qr.hpp:249
pa::LimitedMemoryQR
Incremental QR factorization using modified Gram-Schmidt with reorthogonalization.
Definition: limited-memory-qr.hpp:14
pa::LimitedMemoryQR::q_idx
size_t q_idx
Number of columns of Q being stored.
Definition: limited-memory-qr.hpp:252
pa::LimitedMemoryQR::LimitedMemoryQR
LimitedMemoryQR(size_t n, size_t m)
Definition: limited-memory-qr.hpp:25
pa::LimitedMemoryQR::r_pred
size_t r_pred(size_t i) const
Get the previous index in the circular storage for R.
Definition: limited-memory-qr.hpp:261
pa::LimitedMemoryQR::get_raw_R
const mat & get_raw_R() const
Get the full, raw storage for the upper triangular matrix R.
Definition: limited-memory-qr.hpp:171
pa::LimitedMemoryQR::remove_column
void remove_column()
Remove the leftmost column.
Definition: limited-memory-qr.hpp:75
pa::LimitedMemoryQR::num_columns
size_t num_columns() const
Get the number of columns that are currently stored.
Definition: limited-memory-qr.hpp:227
pa::vec
realvec vec
Default type for vectors.
Definition: vec.hpp:14
pa::LimitedMemoryQR::ring_reverse_iter
ReverseCircularRange< size_t > ring_reverse_iter() const
Get reverse iterators in the circular buffer.
Definition: limited-memory-qr.hpp:244
pa::LimitedMemoryQR::get_raw_Q
const mat & get_raw_Q() const
Get the full, raw storage for the orthogonal matrix Q.
Definition: limited-memory-qr.hpp:166
pa::LimitedMemoryQR::get_reorth_count
size_t get_reorth_count() const
Get the number of MGS reorthogonalizations.
Definition: limited-memory-qr.hpp:206
pa::LimitedMemoryQR::get_Q
mat get_Q() const
Get the matrix Q such that Q times R is the original matrix.
Definition: limited-memory-qr.hpp:197
pa::LimitedMemoryQR::solve
solve_ret_t< Derived > solve(const Eigen::DenseBase< Derived > &B)
Solve the least squares problem AX = B.
Definition: limited-memory-qr.hpp:159
pa::CircularRange
Definition: ringbuffer.hpp:146
panocpy.test.v
v
Definition: test.py:42
pa
Definition: alm.hpp:10
pa::LimitedMemoryQR::R
mat R
Storage for upper triangular factor R.
Definition: limited-memory-qr.hpp:250
panocpy.test.x
x
Definition: test.py:40
pa::LimitedMemoryQR::get_full_R
mat get_full_R() const
Get the full storage for the upper triangular matrix R but with the columns in the correct order.
Definition: limited-memory-qr.hpp:176
bicycle-obstacle-avoidance-mpc.c
c
Definition: bicycle-obstacle-avoidance-mpc.py:124
pa::ReverseCircularRange
Definition: ringbuffer.hpp:181
pa::LimitedMemoryQR::solve_ret_t
std::conditional_t< Eigen::internal::traits< Derived >::ColsAtCompileTime==1, vec, mat > solve_ret_t
Definition: limited-memory-qr.hpp:155
pa::LimitedMemoryQR::ring_head
size_t ring_head() const
Get the head index of the circular buffer (points to the oldest element).
Definition: limited-memory-qr.hpp:230
pa::LimitedMemoryQR::clear_reorth_count
void clear_reorth_count()
Reset the number of MGS reorthogonalizations.
Definition: limited-memory-qr.hpp:208
pa::LimitedMemoryQR::resize
void resize(size_t n, size_t m)
Re-allocate storage for a problem with a different size.
Definition: limited-memory-qr.hpp:220
pa::LimitedMemoryQR::get_R
mat get_R() const
Get the matrix R such that Q times R is the original matrix.
Definition: limited-memory-qr.hpp:190
pa::LimitedMemoryQR::r_succ
size_t r_succ(size_t i) const
Get the next index in the circular storage for R.
Definition: limited-memory-qr.hpp:259
pa::LimitedMemoryQR::n
size_t n() const
Definition: limited-memory-qr.hpp:27
pa::LimitedMemoryQR::r_idx_end
size_t r_idx_end
Index of the one-past-last column of R.
Definition: limited-memory-qr.hpp:254
main.b
b
Definition: main.py:11
pa::LimitedMemoryQR::add_column
void add_column(const VecV &v)
Add the given column to the right.
Definition: limited-memory-qr.hpp:32
pa::LimitedMemoryQR::reset
void reset()
Reset all indices, clearing the Q and R matrices.
Definition: limited-memory-qr.hpp:211
pa::LimitedMemoryQR::ring_next
size_t ring_next(size_t i) const
Get the next index in the circular buffer.
Definition: limited-memory-qr.hpp:235
pa::LimitedMemoryQR::ring_iter
CircularRange< size_t > ring_iter() const
Get iterators in the circular buffer.
Definition: limited-memory-qr.hpp:240
pa::LimitedMemoryQR::ring_tail
size_t ring_tail() const
Get the tail index of the circular buffer (points to one past the most recent element).
Definition: limited-memory-qr.hpp:233
pa::mat
realmat mat
Default type for matrices.
Definition: vec.hpp:20
ringbuffer.hpp
vec.hpp
rosenbrock.X
X
Definition: rosenbrock.py:23
pa::LimitedMemoryQR::solve_col
void solve_col(const VecB &b, VecX &x) const
Solve the least squares problem Ax = b.
Definition: limited-memory-qr.hpp:114
pa::LimitedMemoryQR::r_idx_start
size_t r_idx_start
Index of the first column of R.
Definition: limited-memory-qr.hpp:253
pa::LimitedMemoryQR::scale_R
void scale_R(real_t scal)
Multiply the matrix R by a scalar.
Definition: limited-memory-qr.hpp:200
pa::real_t
double real_t
Default floating point type.
Definition: vec.hpp:8
pa::LimitedMemoryQR::solve
void solve(const MatB &B, MatX &X) const
Solve the least squares problem AX = B.
Definition: limited-memory-qr.hpp:141
pa::LimitedMemoryQR::LimitedMemoryQR
LimitedMemoryQR()=default
pa::LimitedMemoryQR::reorth_count
size_t reorth_count
Number of MGS reorthogonalizations.
Definition: limited-memory-qr.hpp:256
pa::LimitedMemoryQR::ring_prev
size_t ring_prev(size_t i) const
Get the previous index in the circular buffer.
Definition: limited-memory-qr.hpp:237