#include <alpaqa/include/alpaqa/accelerators/internal/limited-memory-qr.hpp>
Incremental QR factorization using modified Gram-Schmidt with reorthogonalization.
Computes A = QR while allowing efficient removal of the first column of A or adding new columns at the end of A.
Definition at line 14 of file limited-memory-qr.hpp.
Public Types | |
template<class Derived > | |
using | solve_ret_t = std::conditional_t< Eigen::internal::traits< Derived >::ColsAtCompileTime==1, vec, mat > |
Public Member Functions | |
LimitedMemoryQR ()=default | |
LimitedMemoryQR (length_t n, length_t m) | |
length_t | n () const |
length_t | m () const |
length_t | size () const |
length_t | history () const |
template<class VecV > | |
void | add_column (const VecV &v) |
Add the given column to the right. | |
void | remove_column () |
Remove the leftmost column. | |
template<class VecB , class VecX > | |
void | solve_col (const VecB &b, VecX &x, real_t tol=0) const |
Solve the least squares problem Ax = b. | |
template<class MatB , class MatX > | |
void | solve (const MatB &B, MatX &X, real_t tol=0) const |
Solve the least squares problem AX = B. | |
template<class Derived > | |
solve_ret_t< Derived > | solve (const Eigen::DenseBase< Derived > &B) |
Solve the least squares problem AX = B. | |
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. | |
mat | get_full_R () const |
Get the full storage for the upper triangular matrix R but with the columns in the correct order. | |
mat | get_R () const |
Get the matrix R such that Q times R is the original matrix. | |
mat | get_Q () const |
Get the matrix Q such that Q times R is the original matrix. | |
void | scale_R (real_t scal) |
Multiply the matrix R by a scalar. | |
unsigned long | get_reorth_count () const |
Get the number of MGS reorthogonalizations. | |
void | clear_reorth_count () |
Reset the number of MGS reorthogonalizations. | |
real_t | get_min_eig () const |
Get the minimum eigenvalue of R. | |
real_t | get_max_eig () const |
Get the maximum eigenvalue of R. | |
void | reset () |
Reset all indices, clearing the Q and R matrices. | |
void | resize (length_t n, length_t m) |
Re-allocate storage for a problem with a different size. | |
length_t | num_columns () const |
Get the number of columns that are currently stored. | |
index_t | ring_head () const |
Get the head index of the circular buffer (points to the oldest element). | |
index_t | ring_tail () const |
Get the tail index of the circular buffer (points to one past the most recent element). | |
index_t | ring_next (index_t i) const |
Get the next index in the circular buffer. | |
index_t | ring_prev (index_t i) const |
Get the previous index in the circular buffer. | |
length_t | current_history () const |
Get the number of columns currently stored in the buffer. | |
CircularRange< index_t > | ring_iter () const |
Get iterators in the circular buffer. | |
ReverseCircularRange< index_t > | ring_reverse_iter () const |
Get reverse iterators in the circular buffer. | |
Private Member Functions | |
index_t | r_succ (index_t i) const |
Get the next index in the circular storage for R. | |
index_t | r_pred (index_t i) const |
Get the previous index in the circular storage for R. | |
Private Attributes | |
mat | Q |
Storage for orthogonal factor Q. | |
mat | R |
Storage for upper triangular factor R. | |
index_t | q_idx = 0 |
Number of columns of Q being stored. | |
index_t | r_idx_start = 0 |
Index of the first column of R. | |
index_t | r_idx_end = 0 |
Index of the one-past-last column of R. | |
unsigned long | reorth_count = 0 |
Number of MGS reorthogonalizations. | |
real_t | min_eig = +inf<config_t> |
Minimum eigenvalue of R. | |
real_t | max_eig = -inf<config_t> |
Maximum eigenvalue of R. | |
using solve_ret_t = std::conditional_t< Eigen::internal::traits<Derived>::ColsAtCompileTime == 1, vec, mat> |
Definition at line 171 of file limited-memory-qr.hpp.
|
default |
|
inline |
n | The size of the vectors, the number of rows of A. |
m | The maximum number of columns of A. |
The maximum dimensions of Q are n×m and the maximum dimensions of R are m×m.
Definition at line 26 of file limited-memory-qr.hpp.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Add the given column to the right.
Definition at line 35 of file limited-memory-qr.hpp.
|
inline |
Remove the leftmost column.
Definition at line 81 of file limited-memory-qr.hpp.
|
inline |
Solve the least squares problem Ax = b.
Do not divide by elements that are smaller in absolute value than tol
.
Definition at line 125 of file limited-memory-qr.hpp.
|
inline |
Solve the least squares problem AX = B.
Do not divide by elements that are smaller in absolute value than tol
.
Definition at line 158 of file limited-memory-qr.hpp.
|
inline |
Solve the least squares problem AX = B.
Definition at line 176 of file limited-memory-qr.hpp.
|
inline |
Get the full, raw storage for the orthogonal matrix Q.
Definition at line 183 of file limited-memory-qr.hpp.
|
inline |
Get the full, raw storage for the upper triangular matrix R.
The columns of this matrix are permuted because it's stored as a circular buffer for efficiently appending columns to the end and popping columns from the front.
Definition at line 188 of file limited-memory-qr.hpp.
|
inline |
Get the full storage for the upper triangular matrix R but with the columns in the correct order.
Definition at line 193 of file limited-memory-qr.hpp.
|
inline |
Get the matrix R such that Q times R is the original matrix.
Definition at line 207 of file limited-memory-qr.hpp.
|
inline |
Get the matrix Q such that Q times R is the original matrix.
Definition at line 214 of file limited-memory-qr.hpp.
|
inline |
Multiply the matrix R by a scalar.
Definition at line 217 of file limited-memory-qr.hpp.
|
inline |
Get the number of MGS reorthogonalizations.
Definition at line 225 of file limited-memory-qr.hpp.
|
inline |
Reset the number of MGS reorthogonalizations.
Definition at line 227 of file limited-memory-qr.hpp.
|
inline |
Get the minimum eigenvalue of R.
Definition at line 230 of file limited-memory-qr.hpp.
|
inline |
Get the maximum eigenvalue of R.
Definition at line 232 of file limited-memory-qr.hpp.
|
inline |
Reset all indices, clearing the Q and R matrices.
Definition at line 235 of file limited-memory-qr.hpp.
Re-allocate storage for a problem with a different size.
Causes a reset.
Definition at line 246 of file limited-memory-qr.hpp.
|
inline |
Get the number of columns that are currently stored.
Definition at line 253 of file limited-memory-qr.hpp.
|
inline |
Get the head index of the circular buffer (points to the oldest element).
Definition at line 256 of file limited-memory-qr.hpp.
|
inline |
Get the tail index of the circular buffer (points to one past the most recent element).
Definition at line 259 of file limited-memory-qr.hpp.
Get the next index in the circular buffer.
Definition at line 261 of file limited-memory-qr.hpp.
Get the previous index in the circular buffer.
Definition at line 263 of file limited-memory-qr.hpp.
|
inline |
Get the number of columns currently stored in the buffer.
Definition at line 265 of file limited-memory-qr.hpp.
|
inline |
Get iterators in the circular buffer.
Definition at line 268 of file limited-memory-qr.hpp.
|
inline |
Get reverse iterators in the circular buffer.
Definition at line 272 of file limited-memory-qr.hpp.
Get the next index in the circular storage for R.
Definition at line 290 of file limited-memory-qr.hpp.
Get the previous index in the circular storage for R.
Definition at line 292 of file limited-memory-qr.hpp.
|
private |
Storage for orthogonal factor Q.
Definition at line 277 of file limited-memory-qr.hpp.
|
private |
Storage for upper triangular factor R.
Definition at line 278 of file limited-memory-qr.hpp.
|
private |
Number of columns of Q being stored.
Definition at line 280 of file limited-memory-qr.hpp.
|
private |
Index of the first column of R.
Definition at line 281 of file limited-memory-qr.hpp.
|
private |
Index of the one-past-last column of R.
Definition at line 282 of file limited-memory-qr.hpp.
|
private |
Number of MGS reorthogonalizations.
Definition at line 284 of file limited-memory-qr.hpp.
Minimum eigenvalue of R.
Definition at line 286 of file limited-memory-qr.hpp.
Maximum eigenvalue of R.
Definition at line 287 of file limited-memory-qr.hpp.