Reference documentation for deal.II version 9.5.0
\(\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 | Private Attributes | List of all members
Householder< number > Class Template Reference

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

Inheritance diagram for Householder< number >:
[legend]

Public Types

using size_type = types::global_dof_index
 

Public Member Functions

 Householder ()=default
 
template<typename number2 >
 Householder (const FullMatrix< number2 > &A)
 
template<typename number2 >
void initialize (const FullMatrix< number2 > &A)
 
template<typename number2 >
double least_squares (Vector< number2 > &dst, const Vector< number2 > &src) const
 
template<typename number2 >
double least_squares (BlockVector< number2 > &dst, const BlockVector< number2 > &src) const
 
template<class VectorType >
void vmult (VectorType &dst, const VectorType &src) const
 
template<class VectorType >
void Tvmult (VectorType &dst, const VectorType &src) const
 

Private Attributes

std::vector< number > diagonal
 
FullMatrix< double > storage
 

Detailed Description

template<typename number>
class Householder< number >

QR-decomposition of a full matrix.

This class computes the QR-decomposition of given matrix by the Householder algorithm. Then, the function least_squares() can be used to compute the vector \(x\) minimizing \(\|Ax-b\|\) for a given vector \(b\). The QR decomposition of \(A\) is useful for this purpose because the minimizer is given by the equation \(x=(A^TA)^{-1}A^Tb=(R^TQ^TQR)^{-1}R^TQ^Tb\) which is easy to compute because \(Q\) is an orthogonal matrix, and consequently \(Q^TQ=I\). Thus, \(x=(R^TR)^{-1}R^TQ^Tb=R^{-1}R^{-T}R^TQ^Tb=R^{-1}Q^Tb\). Furthermore, \(R\) is triangular, so applying \(R^{-1}\) to a vector only involves a backward or forward solve.

Implementation details

The class does not in fact store the \(Q\) and \(R\) factors explicitly as matrices. It does store \(R\), but the \(Q\) factor is stored as the product of Householder reflections of the form \(Q_i = I-v_i v_i^T\) where the vectors \(v_i\) are so that they can be stored in the lower-triangular part of an underlying matrix object, whereas \(R\) is stored in the upper triangular part.

The \(v_i\) vectors and the \(R\) matrix now are in conflict because they both want to use the diagonal entry of the matrix, but we can only store one in these positions, of course. Consequently, the entries \((v_i)_i\) are stored separately in the diagonal member variable.

Note
Instantiations for this template are provided for <float> and <double>; others can be generated in application programs (see the section on Template instantiations in the manual).

Definition at line 80 of file householder.h.

Member Typedef Documentation

◆ size_type

template<typename number >
using Householder< number >::size_type = types::global_dof_index

Declare type of container size type.

Definition at line 86 of file householder.h.

Constructor & Destructor Documentation

◆ Householder() [1/2]

template<typename number >
Householder< number >::Householder ( )
default

Create an empty object.

◆ Householder() [2/2]

template<typename number >
template<typename number2 >
Householder< number >::Householder ( const FullMatrix< number2 > &  A)

Create an object holding the QR-decomposition of the matrix \(A\).

Member Function Documentation

◆ initialize()

template<typename number >
template<typename number2 >
void Householder< number >::initialize ( const FullMatrix< number2 > &  A)

Compute the QR-decomposition of the given matrix \(A\).

This overwrites any previously computed QR decomposition.

◆ least_squares() [1/2]

template<typename number >
template<typename number2 >
double Householder< number >::least_squares ( Vector< number2 > &  dst,
const Vector< number2 > &  src 
) const

Solve the least-squares problem for the right hand side src. The returned scalar value is the Euclidean norm of the approximation error.

  • dst contains the solution of the least squares problem on return.
  • src contains the right hand side b of the least squares problem. It will be changed during the algorithm and is unusable on return.

◆ least_squares() [2/2]

template<typename number >
template<typename number2 >
double Householder< number >::least_squares ( BlockVector< number2 > &  dst,
const BlockVector< number2 > &  src 
) const

This function does the same as the previous one, but for BlockVectors.

◆ vmult()

template<typename number >
template<class VectorType >
void Householder< number >::vmult ( VectorType &  dst,
const VectorType &  src 
) const

A wrapper to least_squares(), implementing the standard MatrixType interface.

◆ Tvmult()

template<typename number >
template<class VectorType >
void Householder< number >::Tvmult ( VectorType &  dst,
const VectorType &  src 
) const

A wrapper to least_squares() that implements multiplication with the transpose matrix.

Member Data Documentation

◆ diagonal

template<typename number >
std::vector<number> Householder< number >::diagonal
private

Storage for the diagonal elements of the orthogonal transformation. See the class documentation for more information.

Definition at line 152 of file householder.h.

◆ storage

template<typename number >
FullMatrix<double> Householder< number >::storage
private

Storage that is internally used for the Householder transformation.

Definition at line 157 of file householder.h.


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