27 size_t n()
const {
return Q.rows(); }
28 size_t m()
const {
return Q.cols(); }
40 for (
size_t i = 0; i <
q_idx; ++i) {
53 while (norm_q < η * norm_v) {
55 for (
size_t i = 0; i <
q_idx; ++i) {
81 Eigen::JacobiRotation<real_t> G;
84 while (r <
q_idx - 1) {
87 G.makeGivens(
R(r,
c),
R(r + 1,
c), &
R(r,
c));
98 R.col(cc).applyOnTheLeft(r, r + 1, G.adjoint());
100 Q.block(0, 0,
Q.rows(),
q_idx).applyOnTheRight(r, r + 1, G);
113 template <
class VecB,
class VecX>
122 for (
auto it_d = rev_bgn; it_d != rev_end; ++it_d) {
124 auto [rR, cR] = *it_d;
127 x(rR) =
Q.col(rR).transpose() *
b;
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);
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());
146 for (Eigen::Index cB = 0; cB < B.cols(); ++cB) {
153 template <
class Derived>
155 Eigen::internal::traits<Derived>::ColsAtCompileTime == 1,
vec,
mat>;
158 template <
class Derived>
182 Eigen::PermutationMatrix<Eigen::Dynamic> P(
R.cols());
184 std::rotate(P.indices().data(), P.indices().data() +
r_idx_start,
185 P.indices().data() + P.size());
193 .triangularView<Eigen::Upper>();
202 R.col(r_idx).topRows(i + 1) *= scal;
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; }
Incremental QR factorization using modified Gram-Schmidt with reorthogonalization.
mat R
Storage for upper triangular factor R.
mat Q
Storage for orthogonal factor Q.
size_t ring_prev(size_t i) const
Get the previous index in the circular buffer.
size_t ring_tail() const
Get the tail index of the circular buffer (points to one past the most recent element).
size_t get_reorth_count() const
Get the number of MGS reorthogonalizations.
mat get_Q() const
Get the matrix Q such that Q times R is the original matrix.
void remove_column()
Remove the leftmost column.
mat get_R() const
Get the matrix R such that Q times R is the original matrix.
void add_column(const VecV &v)
Add the given column to the right.
solve_ret_t< Derived > solve(const Eigen::DenseBase< Derived > &B)
Solve the least squares problem AX = B.
size_t r_pred(size_t i) const
Get the previous index in the circular storage for R.
size_t reorth_count
Number of MGS reorthogonalizations.
void solve_col(const VecB &b, VecX &x) const
Solve the least squares problem Ax = b.
void clear_reorth_count()
Reset the number of MGS reorthogonalizations.
size_t q_idx
Number of columns of Q being stored.
void resize(size_t n, size_t m)
Re-allocate storage for a problem with a different size.
size_t num_columns() const
Get the number of columns that are currently stored.
void solve(const MatB &B, MatX &X) const
Solve the least squares problem AX = B.
ReverseCircularRange< size_t > ring_reverse_iter() const
Get reverse iterators in the circular buffer.
LimitedMemoryQR()=default
std::conditional_t< Eigen::internal::traits< Derived >::ColsAtCompileTime==1, vec, mat > solve_ret_t
size_t r_succ(size_t i) const
Get the next index in the circular storage for R.
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 reset()
Reset all indices, clearing the Q and R matrices.
LimitedMemoryQR(size_t n, size_t m)
void scale_R(real_t scal)
Multiply the matrix R by a scalar.
size_t ring_head() const
Get the head index of the circular buffer (points to the oldest element).
CircularRange< size_t > ring_iter() const
Get 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.
size_t r_idx_start
Index of the first column of R.
size_t ring_next(size_t i) const
Get the next index in the circular buffer.
size_t r_idx_end
Index of the one-past-last column of R.
realmat mat
Default type for matrices.
realvec vec
Default type for vectors.
double real_t
Default floating point type.