Reference documentation for deal.II version GIT relicensing-478-g3275795167 2024-04-24 07:10:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
QR< VectorType > Class Template Reference

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

Inheritance diagram for QR< VectorType >:
Inheritance graph
[legend]

Public Types

using Number = typename VectorType::value_type
 

Public Member Functions

 QR ()
 
virtual ~QR ()=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
 
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)
 

Protected Member Functions

void multiply_with_cols (VectorType &y, const Vector< Number > &x) const
 
void multiply_with_colsT (Vector< Number > &y, const VectorType &x) const
 

Protected Attributes

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
 

Private Member Functions

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

Private Attributes

VectorType tmp
 

Detailed Description

template<typename VectorType>
class QR< VectorType >

A class to compute and store the QR factorization of a matrix represented by a set of column vectors.

The class is design to update a given (possibly empty) QR factorization of a matrix \(A\) (constructed incrementally by providing its columns) due to the addition of a new column vector to \(A\). This is equivalent to constructing an orthonormal basis by the Gram-Schmidt procedure. The class also provides update functionality when the first 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.

See sections 6.5.2-6.5.3 on pp. 335-338 in

@Book{Golub2013,
title = {Matrix computations},
publisher = {Johns Hopkins University Press},
year = {2013},
author = {Golub, Gene H and Van Loan, Charles F},
edition = {4},
}

as well as

@article{Daniel1976,
author = {Daniel, James W and Gragg, Walter Bill and Kaufman, Linda and Stewart, Gilbert W},
title = {{Reorthogonalization and stable algorithms for updating the Gram-Schmidt QR factorization}},
journal = {Mathematics of Computation},
year = {1976},
volume = {30},
number = {136},
pages = {772--795},
}
@Article{Reichel1990,
author = {Reichel, L. and Gragg, W. B.},
title = {{Algorithm 686: FORTRAN Subroutines for Updating the QR Decomposition}},
journal = {ACM Trans. Math. Softw.},
year = {1990},
volume = {16},
number = {4},
pages = {369--377},
month = dec,
issn = {0098-3500},
acmid = {98291},
address = {New York, NY, USA},
doi = {10.1145/98267.98291},
issue_date = {Dec. 1990},
numpages = {9},
publisher = {ACM},
url = {http://doi.acm.org/10.1145/98267.98291},
}

Definition at line 239 of file qr.h.

Member Typedef Documentation

◆ Number

template<typename VectorType >
using QR< VectorType >::Number = typename VectorType::value_type

Number type for R matrix.

Definition at line 245 of file qr.h.

Constructor & Destructor Documentation

◆ QR()

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

Default constructor.

◆ ~QR()

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

Destructor.

Member Function Documentation

◆ append_column()

template<typename VectorType >
virtual bool QR< 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.

Note
Currently this function always returns true.

Implements BaseQR< VectorType >.

◆ remove_column()

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

Remove first column and update QR factorization.

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

The standard approach is to partition \(R\) as

\[ R = \begin{bmatrix} r_{11} & w^T \\ 0 & R_{33} \end{bmatrix} \]

It then follows that

\[ Q^T \tilde A = \begin{bmatrix} 0 & w^T \\ 0 & R_{33} \end{bmatrix} \]

is upper Hessenberg where unwanted sub-diagonal elements can be zeroed by a sequence of Givens rotations.

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 QR< 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 QR< 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 QR< 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 QR< 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 >.

◆ apply_givens_rotation()

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

Apply givens rotation in the (i,j)-plane to Q and R so that R(k,k) is zeroed.

See Chapter 5.1.9 of Golub 2013, Matrix computations.

◆ size()

template<typename VectorType >
unsigned int BaseQR< VectorType >::size ( ) const
inherited

Return size of the subspace.

◆ get_R()

template<typename VectorType >
const LAPACKFullMatrix< Number > & BaseQR< VectorType >::get_R ( ) const
inherited

Return the current upper triangular matrix R.

◆ solve()

template<typename VectorType >
void BaseQR< VectorType >::solve ( Vector< Number > &  x,
const Vector< Number > &  y,
const bool  transpose = false 
) const
inherited

Solve \(Rx=y\). Vectors x and y should be consistent with the current size of the subspace. If transpose is true, \(R^Tx=y\) is solved instead.

◆ connect_givens_slot()

template<typename VectorType >
boost::signals2::connection BaseQR< VectorType >::connect_givens_slot ( const std::function< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &csr)> &  slot)
inherited

Connect a slot to retrieve a notification when the Givens rotations are performed.

The function takes two indices, i and j, describing the plane of rotation, and a triplet of numbers csr (cosine, sine and radius, see Utilities::LinearAlgebra::givens_rotation()) which represents the rotation matrix.

◆ multiply_with_cols()

template<typename VectorType >
void BaseQR< VectorType >::multiply_with_cols ( VectorType &  y,
const Vector< Number > &  x 
) const
protectedinherited

Compute \(y=Hx\) where \(H\) is the matrix formed by the column vectors stored by this object.

◆ multiply_with_colsT()

template<typename VectorType >
void BaseQR< VectorType >::multiply_with_colsT ( Vector< Number > &  y,
const VectorType &  x 
) const
protectedinherited

Multiply with transpose columns stored in the object.

Member Data Documentation

◆ tmp

template<typename VectorType >
VectorType QR< VectorType >::tmp
private

Temporary vector needed to do Givens rotation of Q.

Definition at line 325 of file qr.h.

◆ columns

template<typename VectorType >
std::vector<std::unique_ptr<VectorType> > BaseQR< VectorType >::columns
protectedinherited

A vector of unique pointers to store columns.

Definition at line 159 of file qr.h.

◆ R

template<typename VectorType >
LAPACKFullMatrix<Number> BaseQR< VectorType >::R
protectedinherited

Matrix to store R.

Definition at line 164 of file qr.h.

◆ current_size

template<typename VectorType >
unsigned int BaseQR< VectorType >::current_size
protectedinherited

Current size (number of columns in Q).

Definition at line 169 of file qr.h.

◆ givens_signal

template<typename VectorType >
boost::signals2::signal<void(const unsigned int i, const unsigned int j, const std::array<Number, 3> &)> BaseQR< VectorType >::givens_signal
protectedinherited

Signal used to retrieve a notification when Givens rotations are performed in the (i,j)-plane.

Definition at line 178 of file qr.h.


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