quala 0.0.1a1
Quasi-Newton and other accelerators
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
LimitedMemoryQR Class Reference

#include <quala/detail/limited-memory-qr.hpp>

Detailed Description

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.

+ Collaboration diagram for LimitedMemoryQR:

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. More...
 
void remove_column ()
 Remove the leftmost column. More...
 
template<class VecB , class VecX >
void solve_col (const VecB &b, VecX &x) const
 Solve the least squares problem Ax = b. More...
 
template<class MatB , class MatX >
void solve (const MatB &B, MatX &X) const
 Solve the least squares problem AX = B. More...
 
template<class Derived >
solve_ret_t< Derived > solve (const Eigen::DenseBase< Derived > &B)
 Solve the least squares problem AX = B. More...
 
const matget_raw_Q () const
 Get the full, raw storage for the orthogonal matrix Q. More...
 
const matget_raw_R () const
 Get the full, raw storage for the upper triangular matrix R. More...
 
mat get_full_R () const
 Get the full storage for the upper triangular matrix R but with the columns in the correct order. More...
 
mat get_R () const
 Get the matrix R such that Q times R is the original matrix. More...
 
mat get_Q () const
 Get the matrix Q such that Q times R is the original matrix. More...
 
void scale_R (real_t scal)
 Multiply the matrix R by a scalar. More...
 
unsigned long get_reorth_count () const
 Get the number of MGS reorthogonalizations. More...
 
void clear_reorth_count ()
 Reset the number of MGS reorthogonalizations. More...
 
void reset ()
 Reset all indices, clearing the Q and R matrices. More...
 
void resize (length_t n, length_t m)
 Re-allocate storage for a problem with a different size. More...
 
length_t num_columns () const
 Get the number of columns that are currently stored. More...
 
index_t ring_head () const
 Get the head index of the circular buffer (points to the oldest element). More...
 
index_t ring_tail () const
 Get the tail index of the circular buffer (points to one past the most recent element). More...
 
index_t ring_next (index_t i) const
 Get the next index in the circular buffer. More...
 
index_t ring_prev (index_t i) const
 Get the previous index in the circular buffer. More...
 
length_t current_history () const
 Get the number of columns currently stored in the buffer. More...
 
CircularRange< index_tring_iter () const
 Get iterators in the circular buffer. More...
 
ReverseCircularRange< index_tring_reverse_iter () const
 Get reverse iterators in the circular buffer. More...
 

Private Member Functions

index_t r_succ (index_t i) const
 Get the next index in the circular storage for R. More...
 
index_t r_pred (index_t i) const
 Get the previous index in the circular storage for R. More...
 

Private Attributes

mat Q
 Storage for orthogonal factor Q. More...
 
mat R
 Storage for upper triangular factor R. More...
 
index_t q_idx = 0
 Number of columns of Q being stored. More...
 
index_t r_idx_start = 0
 Index of the first column of R. More...
 
index_t r_idx_end = 0
 Index of the one-past-last column of R. More...
 
unsigned long reorth_count = 0
 Number of MGS reorthogonalizations. More...
 

Member Typedef Documentation

◆ solve_ret_t

using solve_ret_t = std::conditional_t< Eigen::internal::traits<Derived>::ColsAtCompileTime == 1, vec, mat>

Definition at line 156 of file limited-memory-qr.hpp.

Constructor & Destructor Documentation

◆ LimitedMemoryQR() [1/2]

LimitedMemoryQR ( )
default

◆ LimitedMemoryQR() [2/2]

LimitedMemoryQR ( length_t  n,
length_t  m 
)
inline
Parameters
nThe size of the vectors, the number of rows of A.
mThe 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 25 of file limited-memory-qr.hpp.

Member Function Documentation

◆ n()

length_t n ( ) const
inline

Definition at line 27 of file limited-memory-qr.hpp.

+ Here is the caller graph for this function:

◆ m()

length_t m ( ) const
inline

Definition at line 28 of file limited-memory-qr.hpp.

+ Here is the caller graph for this function:

◆ size()

length_t size ( ) const
inline

Definition at line 29 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ history()

length_t history ( ) const
inline

Definition at line 30 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ add_column()

void add_column ( const VecV &  v)
inline

Add the given column to the right.

Definition at line 34 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ remove_column()

void remove_column ( )
inline

Remove the leftmost column.

Definition at line 77 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve_col()

void solve_col ( const VecB &  b,
VecX &  x 
) const
inline

Solve the least squares problem Ax = b.

Definition at line 116 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve() [1/2]

void solve ( const MatB &  B,
MatX &  X 
) const
inline

Solve the least squares problem AX = B.

Definition at line 143 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve() [2/2]

solve_ret_t< Derived > solve ( const Eigen::DenseBase< Derived > &  B)
inline

Solve the least squares problem AX = B.

Definition at line 161 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ get_raw_Q()

const mat & get_raw_Q ( ) const
inline

Get the full, raw storage for the orthogonal matrix Q.

Definition at line 168 of file limited-memory-qr.hpp.

◆ get_raw_R()

const mat & get_raw_R ( ) const
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 173 of file limited-memory-qr.hpp.

◆ get_full_R()

mat get_full_R ( ) const
inline

Get the full storage for the upper triangular matrix R but with the columns in the correct order.

Note
Meant for tests only, creates a permuted copy.

Definition at line 178 of file limited-memory-qr.hpp.

+ Here is the caller graph for this function:

◆ get_R()

mat get_R ( ) const
inline

Get the matrix R such that Q times R is the original matrix.

Note
Meant for tests only, creates a permuted copy.

Definition at line 192 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ get_Q()

mat get_Q ( ) const
inline

Get the matrix Q such that Q times R is the original matrix.

Note
Meant for tests only, creates a copy.

Definition at line 199 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ scale_R()

void scale_R ( real_t  scal)
inline

Multiply the matrix R by a scalar.

Definition at line 202 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ get_reorth_count()

unsigned long get_reorth_count ( ) const
inline

Get the number of MGS reorthogonalizations.

Definition at line 208 of file limited-memory-qr.hpp.

◆ clear_reorth_count()

void clear_reorth_count ( )
inline

Reset the number of MGS reorthogonalizations.

Definition at line 210 of file limited-memory-qr.hpp.

◆ reset()

void reset ( )
inline

Reset all indices, clearing the Q and R matrices.

Definition at line 213 of file limited-memory-qr.hpp.

+ Here is the caller graph for this function:

◆ resize()

void resize ( length_t  n,
length_t  m 
)
inline

Re-allocate storage for a problem with a different size.

Causes a reset.

Definition at line 222 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ num_columns()

length_t num_columns ( ) const
inline

Get the number of columns that are currently stored.

Definition at line 229 of file limited-memory-qr.hpp.

+ Here is the caller graph for this function:

◆ ring_head()

index_t ring_head ( ) const
inline

Get the head index of the circular buffer (points to the oldest element).

Definition at line 232 of file limited-memory-qr.hpp.

◆ ring_tail()

index_t ring_tail ( ) const
inline

Get the tail index of the circular buffer (points to one past the most recent element).

Definition at line 235 of file limited-memory-qr.hpp.

+ Here is the caller graph for this function:

◆ ring_next()

index_t ring_next ( index_t  i) const
inline

Get the next index in the circular buffer.

Definition at line 237 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ ring_prev()

index_t ring_prev ( index_t  i) const
inline

Get the previous index in the circular buffer.

Definition at line 239 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:

◆ current_history()

length_t current_history ( ) const
inline

Get the number of columns currently stored in the buffer.

Definition at line 241 of file limited-memory-qr.hpp.

+ Here is the caller graph for this function:

◆ ring_iter()

CircularRange< index_t > ring_iter ( ) const
inline

Get iterators in the circular buffer.

Definition at line 244 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ring_reverse_iter()

ReverseCircularRange< index_t > ring_reverse_iter ( ) const
inline

Get reverse iterators in the circular buffer.

Definition at line 248 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ r_succ()

index_t r_succ ( index_t  i) const
inlineprivate

Get the next index in the circular storage for R.

Definition at line 263 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ r_pred()

index_t r_pred ( index_t  i) const
inlineprivate

Get the previous index in the circular storage for R.

Definition at line 265 of file limited-memory-qr.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Member Data Documentation

◆ Q

mat Q
private

Storage for orthogonal factor Q.

Definition at line 253 of file limited-memory-qr.hpp.

◆ R

mat R
private

Storage for upper triangular factor R.

Definition at line 254 of file limited-memory-qr.hpp.

◆ q_idx

index_t q_idx = 0
private

Number of columns of Q being stored.

Definition at line 256 of file limited-memory-qr.hpp.

◆ r_idx_start

index_t r_idx_start = 0
private

Index of the first column of R.

Definition at line 257 of file limited-memory-qr.hpp.

◆ r_idx_end

index_t r_idx_end = 0
private

Index of the one-past-last column of R.

Definition at line 258 of file limited-memory-qr.hpp.

◆ reorth_count

unsigned long reorth_count = 0
private

Number of MGS reorthogonalizations.

Definition at line 260 of file limited-memory-qr.hpp.


The documentation for this class was generated from the following file: