Public Types | Public Member Functions | Protected Member Functions | List of all members
QRMatrix< MatrixType > Class Template Reference

QRMatrix (i.e. QR decomposition, QR factorisation or orthogonal-triangular decomposition) decomposes a scalar/complex matrix A into the following matrix product: More...

Public Types

enum  outputTypes : uint8_t { FULL_R = 1, FULL_QR = 2, REDUCED_R = 3 }
 
enum  storeMethods : uint8_t { IN_PLACE = 1, OUT_OF_PLACE = 2 }
 
enum  colPivoting : bool { FALSE = false, TRUE = true }
 
typedef MatrixType::cmptType cmptType
 
typedef SquareMatrix< cmptTypeSMatrix
 
typedef RectangularMatrix< cmptTypeRMatrix
 

Public Member Functions

 QRMatrix ()
 
 QRMatrix (const outputTypes outputType, const storeMethods storeMethod=storeMethods::OUT_OF_PLACE, const colPivoting colPivot=colPivoting::FALSE)
 
 QRMatrix (MatrixType &A, const outputTypes outputType, const storeMethods storeMethod=storeMethods::OUT_OF_PLACE, const colPivoting colPivot=colPivoting::FALSE)
 
 QRMatrix (const MatrixType &A, const outputTypes outputType, const colPivoting colPivot=colPivoting::FALSE)
 
const MatrixType & Q () const
 
const MatrixType & R () const
 
const labelListorderP () const
 
SMatrix P () const
 
void decompose (MatrixType &A)
 
void decompose (const MatrixType &A)
 
void solve (List< cmptType > &x, const UList< cmptType > &source) const
 
template<class Addr >
void solve (List< cmptType > &x, const IndirectListBase< cmptType, Addr > &source) const
 
tmp< Field< cmptType > > solve (const UList< cmptType > &source) const
 
template<class Addr >
tmp< Field< cmptType > > solve (const IndirectListBase< cmptType, Addr > &source) const
 
SMatrix inv () const
 
RMatrix backSubstitution (const SMatrix &A, const RMatrix &b)
 
template<class Addr >
Foam::tmp< Foam::Field< typename MatrixType::cmptType > > solve (const IndirectListBase< cmptType, Addr > &source) const
 

Protected Member Functions

RMatrix householderReflector (RMatrix u)
 
void applyLeftReflector (MatrixType &A, const RMatrix &reflector, const label k=0, const label k1=0)
 
void applyRightReflector (MatrixType &A, const RMatrix &reflector, const label k=0)
 

Detailed Description

template<class MatrixType>
class Foam::QRMatrix< MatrixType >

QRMatrix (i.e. QR decomposition, QR factorisation or orthogonal-triangular decomposition) decomposes a scalar/complex matrix A into the following matrix product:

    A = Q*R,

where Q is a unitary similarity matrix, R is an upper triangular matrix.

Reference:

    QR decomposition:
        mathworld.wolfram.com/QRDecomposition.html (Retrieved:17-06-19)
        w.wiki/52X (Retrieved: 17-06-19)

    QR decomposition with Householder reflector (tag:M):
        Monahan, J. F. (2011).
        Numerical methods of statistics.
        Cambridge: Cambridge University Press.
        DOI:10.1017/CBO9780511977176

    QR decomposition with column pivoting (tag:QSB):
        Quintana-Ortí, G., Sun, X., & Bischof, C. H. (1998).
        A BLAS-3 version of the QR factorization with column pivoting.
        SIAM Journal on Scientific Computing, 19(5), 1486-1494.
        DOI:10.1137/S1064827595296732

    Moore-Penrose inverse algorithm (tags:KP; KPP):
        Katsikis, V. N., & Pappas, D. (2008).
        Fast computing of the Moore-Penrose inverse matrix.
        Electronic Journal of Linear Algebra, 17(1), 637-650.
        DOI:10.13001/1081-3810.1287

        Katsikis, V. N., Pappas, D., & Petralias, A. (2011).
        An improved method for the computation of
        the Moore–Penrose inverse matrix.
        Applied Mathematics and Computation, 217(23), 9828-9834.
        DOI:10.1016/j.amc.2011.04.080

    Tolerance for the Moore-Penrose inverse algorithm (tag:TA):
        Toutounian, F., & Ataei, A. (2009).
        A new method for computing Moore–Penrose inverse matrices.
        Journal of Computational and applied Mathematics, 228(1), 412-417.
        DOI:10.1016/j.cam.2008.10.008

    Operands of QR decomposition:
        mathworld.wolfram.com/UnitaryMatrix.html (Retrieved:13-06-19)
        mathworld.wolfram.com/UpperTriangularMatrix.html (Retrieved:16-06-19)
        mathworld.wolfram.com/PermutationMatrix.html (Retrieved:20-06-19)

    Back substitution:
        mathworld.wolfram.com/GaussianElimination.html (Retrieved:15-06-19)
Usage
Input types:
  • A can be a SquareMatrix<Type> or RectangularMatrix<Type>

Output types:

Options for the output forms of QRMatrix (for an (m-by-n) input matrix A with k = min(m, n)):

Options where to store R:

Options for the computation of column pivoting:

Direct solution of linear systems A x = b is possible by solve() alongside the following limitations:

Notes

See also
Test-QRMatrix.C
Source files

Definition at line 148 of file QRMatrix.H.

Member Typedef Documentation

◆ cmptType

typedef MatrixType::cmptType cmptType

Definition at line 153 of file QRMatrix.H.

◆ SMatrix

Definition at line 154 of file QRMatrix.H.

◆ RMatrix

Definition at line 155 of file QRMatrix.H.

Member Enumeration Documentation

◆ outputTypes

enum outputTypes : uint8_t
Enumerator
FULL_R 

computes only R

FULL_QR 

computes both R and Q

REDUCED_R 

computes only reduced R

Definition at line 158 of file QRMatrix.H.

◆ storeMethods

enum storeMethods : uint8_t
Enumerator
IN_PLACE 

replaces input matrix content with R

OUT_OF_PLACE 

creates new object of R

Definition at line 166 of file QRMatrix.H.

◆ colPivoting

enum colPivoting : bool
Enumerator
FALSE 

switches off column pivoting

TRUE 

switches on column pivoting

Definition at line 173 of file QRMatrix.H.

Constructor & Destructor Documentation

◆ QRMatrix() [1/4]

Definition at line 179 of file QRMatrix.C.

◆ QRMatrix() [2/4]

QRMatrix ( const outputTypes  outputType,
const storeMethods  storeMethod = storeMethods::OUT_OF_PLACE,
const colPivoting  colPivot = colPivoting::FALSE 
)

Definition at line 189 of file QRMatrix.C.

◆ QRMatrix() [3/4]

QRMatrix ( MatrixType &  A,
const outputTypes  outputType,
const storeMethods  storeMethod = storeMethods::OUT_OF_PLACE,
const colPivoting  colPivot = colPivoting::FALSE 
)

Definition at line 206 of file QRMatrix.C.

References A.

◆ QRMatrix() [4/4]

QRMatrix ( const MatrixType &  A,
const outputTypes  outputType,
const colPivoting  colPivot = colPivoting::FALSE 
)

Definition at line 226 of file QRMatrix.C.

References A.

Member Function Documentation

◆ householderReflector()

Foam::RectangularMatrix< typename MatrixType::cmptType > householderReflector ( RMatrix  u)
inlineprotected

Definition at line 45 of file QRMatrixI.H.

References Foam::abort(), Matrix::columnNorm(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::mag(), and Matrix::n().

Here is the call graph for this function:

◆ applyLeftReflector()

void applyLeftReflector ( MatrixType &  A,
const RMatrix reflector,
const label  k = 0,
const label  k1 = 0 
)
inlineprotected

Definition at line 90 of file QRMatrixI.H.

References A, Foam::Detail::conj(), k, Matrix::m(), Foam::sum(), and Foam::Zero.

Here is the call graph for this function:

◆ applyRightReflector()

void applyRightReflector ( MatrixType &  A,
const RMatrix reflector,
const label  k = 0 
)
inlineprotected

Definition at line 119 of file QRMatrixI.H.

References A, Foam::Detail::conj(), k, Matrix::m(), Foam::sum(), and Foam::Zero.

Here is the call graph for this function:

◆ Q()

const MatrixType & Q
inline

Definition at line 148 of file QRMatrixI.H.

Referenced by Foam::pinv().

Here is the caller graph for this function:

◆ R()

const MatrixType & R
inline

Definition at line 157 of file QRMatrixI.H.

Referenced by Foam::pinv().

Here is the caller graph for this function:

◆ orderP()

const Foam::labelList & orderP
inline

Definition at line 164 of file QRMatrixI.H.

◆ P()

Foam::SquareMatrix< typename Foam::QRMatrix< MatrixType >::cmptType > P
inline

Definition at line 172 of file QRMatrixI.H.

References forAll, and Foam::Zero.

Referenced by Foam::pinv().

Here is the caller graph for this function:

◆ decompose() [1/2]

void decompose ( MatrixType &  A)

Definition at line 247 of file QRMatrix.C.

References A, Foam::I, Foam::nl, and WarningInFunction.

◆ decompose() [2/2]

void decompose ( const MatrixType &  A)

Definition at line 306 of file QRMatrix.C.

References A, Foam::I, Foam::nl, and WarningInFunction.

◆ solve() [1/5]

void solve ( List< cmptType > &  x,
const UList< cmptType > &  source 
) const

Definition at line 351 of file QRMatrix.C.

References x.

◆ solve() [2/5]

void solve ( List< cmptType > &  x,
const IndirectListBase< cmptType, Addr > &  source 
) const

Definition at line 363 of file QRMatrix.C.

References x.

◆ solve() [3/5]

Foam::tmp< Foam::Field< typename MatrixType::cmptType > > solve ( const UList< cmptType > &  source) const

Definition at line 375 of file QRMatrix.C.

◆ solve() [4/5]

tmp<Field<cmptType> > solve ( const IndirectListBase< cmptType, Addr > &  source) const

◆ inv()

Foam::QRMatrix< MatrixType >::SMatrix inv

Definition at line 405 of file QRMatrix.C.

References Foam::inv(), Matrix::m(), and x.

Here is the call graph for this function:

◆ backSubstitution()

Foam::RectangularMatrix< typename MatrixType::cmptType > backSubstitution ( const SMatrix A,
const RMatrix b 
)

Definition at line 429 of file QRMatrix.C.

References A, Foam::abort(), Foam::constant::atomic::alpha, Foam::diag(), Foam::FatalError, FatalErrorInFunction, k, Matrix::m(), Foam::mag(), Matrix::n(), Foam::nl, Foam::tab, WarningInFunction, and Foam::Zero.

Referenced by Foam::pinv().

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

◆ solve() [5/5]

Foam::tmp<Foam::Field<typename MatrixType::cmptType> > solve ( const IndirectListBase< cmptType, Addr > &  source) const

Definition at line 391 of file QRMatrix.C.


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