13template <Config Conf = DefaultConfig>
56 while (norm_q < η * norm_v) {
87 Eigen::JacobiRotation<real_t> G;
91 while (r <
q_idx - 1) {
94 G.makeGivens(
R(r, c),
R(r + 1, c), &
R(r, c));
105 R.col(cc).applyOnTheLeft(r, r + 1, G.adjoint());
107 Q.block(0, 0,
Q.rows(),
q_idx).applyOnTheRight(r, r + 1, G);
124 template <
class VecB,
class VecX>
133 for (
auto it_d = rev_bgn; it_d != rev_end; ++it_d) {
135 auto [rR, cR] = *it_d;
137 if (std::abs(
R(rR, cR)) < tol) {
143 x(rR) =
Q.col(rR).transpose() * b;
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);
157 template <
class MatB,
class MatX>
159 assert(B.cols() <= X.cols());
160 assert(B.rows() >=
Q.rows());
163 for (Eigen::Index cB = 0; cB < B.cols(); ++cB) {
170 template <
class Derived>
172 Eigen::internal::traits<Derived>::ColsAtCompileTime == 1,
vec,
mat>;
175 template <
class Derived>
199 Eigen::PermutationMatrix<Eigen::Dynamic> P(
R.cols());
201 std::rotate(P.indices().data(), P.indices().data() +
r_idx_start,
202 P.indices().data() + P.size());
210 .template triangularView<Eigen::Upper>();
219 R.col(r_idx).topRows(i + 1) *= scal;
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.
LimitedMemoryQR()=default
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)
typename Conf::real_t real_t
typename Conf::index_t index_t
typename Conf::length_t length_t