Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
ImplicitQR< VectorType > Class Template Reference

#include <deal.II/lac/qr.h>

Inheritance diagram for ImplicitQR< VectorType >:
[legend]

Public Types

typedef VectorType::value_type Number
 

Public Member Functions

 ImplicitQR ()
 
virtual ~ImplicitQR ()=default
 
virtual bool append_column (const VectorType &column)
 
virtual void remove_column (const unsigned int k=0)
 
virtual void multiply_with_Q (VectorType &y, const Vector< Number > &x) const
 
virtual void multiply_with_QT (Vector< Number > &y, const VectorType &x) const
 
virtual void multiply_with_A (VectorType &y, const Vector< Number > &x) const
 
virtual void multiply_with_AT (Vector< Number > &y, const VectorType &x) const
 
boost::signals2::connection connect_append_column_slot (const std::function< bool(const Vector< Number > &u, const Number &rho2, const Number &col_norm_sqr)> &slot)
 
- Public Member Functions inherited from BaseQR< VectorType >
virtual ~BaseQR ()=default
 
unsigned int size () const
 
const LAPACKFullMatrix< Number > & get_R () const
 
void solve (Vector< Number > &x, const Vector< Number > &y, const bool transpose=false) const
 
boost::signals2::connection connect_givens_slot (const std::function< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &csr)> &slot)
 

Private Member Functions

void apply_givens_rotation (const unsigned int i, const unsigned int k)
 

Private Attributes

boost::signals2::signal< bool(const Vector< Number > &u, const Number &rho, const Number &col_norm_sqr)> column_signal
 

Additional Inherited Members

- Protected Member Functions inherited from BaseQR< VectorType >
 BaseQR ()
 
void multiply_with_cols (VectorType &y, const Vector< Number > &x) const
 
void multiply_with_colsT (Vector< Number > &y, const VectorType &x) const
 
- Protected Attributes inherited from BaseQR< VectorType >
std::vector< std::unique_ptr< VectorType > > columns
 
LAPACKFullMatrix< NumberR
 
unsigned int current_size
 
boost::signals2::signal< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &)> givens_signal
 

Detailed Description

template<typename VectorType>
class ImplicitQR< VectorType >

A class to obtain the triangular \(R\) matrix of the \(A=QR\) factorization together with the matrix \(A\) itself. The orthonormal matrix \(Q\) is not stored explicitly, the name of the class. The multiplication with \(Q\) can be represented as \(Q=A R^{-1}\), whereas the multiplication with \(Q^T\) is given by \(Q^T=R^{-T}A^T\).

The class is designed to update a given (possibly empty) QR factorization due to the addition of a new column vector. This is equivalent to constructing an orthonormal basis by the Gram-Schmidt procedure. The class also provides update functionality when the column is removed.

The VectorType template argument may either be a parallel and serial vector, and only need to have basic operations such as additions, scalar product, etc. It also needs to have a copy-constructor.

Definition at line 347 of file qr.h.

Member Typedef Documentation

◆ Number

template<typename VectorType >
typedef VectorType::value_type ImplicitQR< VectorType >::Number

Number type for R matrix.

Definition at line 353 of file qr.h.

Constructor & Destructor Documentation

◆ ImplicitQR()

template<typename VectorType >
ImplicitQR< VectorType >::ImplicitQR ( )

Default constructor.

◆ ~ImplicitQR()

template<typename VectorType >
virtual ImplicitQR< VectorType >::~ImplicitQR ( )
virtualdefault

Destructor.

Member Function Documentation

◆ append_column()

template<typename VectorType >
virtual bool ImplicitQR< VectorType >::append_column ( const VectorType &  column)
virtual

Append column to the QR factorization. Returns true if the result is successful, i.e. the columns are linearly independent. Otherwise the column is rejected and the return value is false.

Implements BaseQR< VectorType >.

◆ remove_column()

template<typename VectorType >
virtual void ImplicitQR< VectorType >::remove_column ( const unsigned int  k = 0)
virtual

Remove column and update QR factorization.

Starting from the given QR decomposition \(QR= A = [a_1\,\dots a_n], \quad a_i \in R^m\) we aim at computing factorization of \(\tilde Q \tilde R= \tilde A = [a_2\,\dots a_n], \quad a_i \in R^m\).

Note that \(\tilde R^T \tilde R = \tilde A^T \tilde A\), where the RHS is included in \(A^T A = R^T R\). Therefore \(\tilde R\) can be obtained by Cholesky decomposition.

Implements BaseQR< VectorType >.

◆ multiply_with_Q()

template<typename VectorType >
virtual void ImplicitQR< VectorType >::multiply_with_Q ( VectorType &  y,
const Vector< Number > &  x 
) const
virtual

Set \(y = Qx\). The size of \(x\) should be consistent with the size of the R matrix.

Implements BaseQR< VectorType >.

◆ multiply_with_QT()

template<typename VectorType >
virtual void ImplicitQR< VectorType >::multiply_with_QT ( Vector< Number > &  y,
const VectorType &  x 
) const
virtual

Set \(y = Q^Tx\). The size of \(x\) should be consistent with the size of column vectors.

Implements BaseQR< VectorType >.

◆ multiply_with_A()

template<typename VectorType >
virtual void ImplicitQR< VectorType >::multiply_with_A ( VectorType &  y,
const Vector< Number > &  x 
) const
virtual

Set \(y = QRx\). The size of \(x\) should be consistent with the size of the R matrix.

Implements BaseQR< VectorType >.

◆ multiply_with_AT()

template<typename VectorType >
virtual void ImplicitQR< VectorType >::multiply_with_AT ( Vector< Number > &  y,
const VectorType &  x 
) const
virtual

Set \(y = R^T Q^Tx\). The size of \(x\) should be consistent with the size of column vectors.

Implements BaseQR< VectorType >.

◆ connect_append_column_slot()

template<typename VectorType >
boost::signals2::connection ImplicitQR< VectorType >::connect_append_column_slot ( const std::function< bool(const Vector< Number > &u, const Number &rho2, const Number &col_norm_sqr)> &  slot)

Connect a slot to implement a custom check of linear dependency during addition of a column.

Here, u is the last column of the to-be R matrix, rho is its diagonal and col_norm_sqr is the square of the \(l2\) norm of the column. The function should return true if the new column is linearly independent.

◆ apply_givens_rotation()

template<typename VectorType >
void ImplicitQR< VectorType >::apply_givens_rotation ( const unsigned int  i,
const unsigned int  k 
)
private

Apply givens rotation in the (i,k)-plane to zero out \(R(k,k)\).

Member Data Documentation

◆ column_signal

template<typename VectorType >
boost::signals2::signal<bool(const Vector<Number> &u, const Number & rho, const Number & col_norm_sqr)> ImplicitQR< VectorType >::column_signal
private

Signal used to decide if the new column is linear dependent.

Here, u is the last column of the to-be R matrix, rho is its diagonal and col_norm_sqr is the square of the \(l2\) norm of the column. The function should return true if the new column is linearly independent.

Definition at line 428 of file qr.h.


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