QRMatrix< MatrixType > Class Template Reference

QRMatrix computes QR decomposition of a given scalar/complex matrix A into the following: More...

Public Types

enum  modes : uint8_t { FULL = 1, ECONOMY = 2 }
 Options for the decomposition mode. More...
 
enum  outputs : uint8_t { ONLY_Q = 1, ONLY_R = 2, BOTH_QR = 3 }
 Options for the output types. More...
 
enum  pivoting : bool { FALSE = false, TRUE = true }
 Options for the column pivoting. More...
 
typedef MatrixType::cmptType cmptType
 
typedef SquareMatrix< cmptTypeSMatrix
 
typedef RectangularMatrix< cmptTypeRMatrix
 

Public Member Functions

 QRMatrix ()=delete
 Construct null. More...
 
 QRMatrix (const modes mode, const outputs output, const bool pivoting, MatrixType &A)
 Construct with a matrix and perform QR decomposition. More...
 
 QRMatrix (const MatrixType &A, const modes mode, const outputs output=outputs::BOTH_QR, const bool pivoting=false)
 Construct with a const matrix and perform QR decomposition. More...
 
const MatrixType & Q () const noexcept
 Return const reference to Q. More...
 
MatrixType & Q () noexcept
 Return reference to Q. More...
 
const MatrixType & R () const noexcept
 Return const reference to R. More...
 
MatrixType & R () noexcept
 Return reference to R. More...
 
const labelListp () const noexcept
 Return const reference to p. More...
 
SMatrix P () const
 Create and return the permutation matrix. More...
 
void solve (List< cmptType > &x, const UList< cmptType > &source) const
 Solve the linear system with the given source and return the solution in the argument x. More...
 
template<class Addr >
void solve (List< cmptType > &x, const IndirectListBase< cmptType, Addr > &source) const
 Solve the linear system with the given source and return the solution in the argument x. More...
 
tmp< Field< cmptType > > solve (const UList< cmptType > &source) const
 Solve the linear system with the given source and return the solution. More...
 
template<class Addr >
tmp< Field< cmptType > > solve (const IndirectListBase< cmptType, Addr > &source) const
 Solve the linear system with the given source and return the solution. More...
 
RMatrix solve (const RMatrix &b)
 Solve a row-echelon-form linear system (Ax = b) starting from the bottom by back substitution. More...
 
SMatrix inv () const
 Return the inverse of (Q*R), solving x = (Q*R).inv()*source. More...
 
template<class Addr >
Foam::tmp< Foam::Field< typename MatrixType::cmptType > > solve (const IndirectListBase< cmptType, Addr > &source) const
 

Detailed Description

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

QRMatrix computes QR decomposition of a given scalar/complex matrix A into the following:

    A = Q R

or in case of QR decomposition with column pivoting:

    A P = Q R

where

$ Q $ = Unitary/orthogonal matrix
$ R $ = Upper triangular matrix
$ P $ = Permutation matrix

References:

        TNT implementation:
            Pozo, R. (1997).
            Template Numerical Toolkit for linear algebra:
            High performance programming with C++
            and the Standard Template Library.
            The International Journal of Supercomputer Applications
            and High Performance Computing, 11(3), 251-263.
            DOI:10.1177/109434209701100307

        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
Usage


Input:

$ A $ = RectangularMatrix<Type> or SquareMatrix<Type>


Options for the decomposition mode:

$ modes::FULL $ = compute full-size decomposition
$ modes::ECONOMY $ = compute economy-size decomposition


Options for the output types:

$ outputs::ONLY\_Q $ = compute only Q
$ outputs::ONLY\_R $ = compute only R
$ outputs::BOTH\_QR $ = compute both Q and R


Options for the column pivoting:

$ pivoting::FALSE $ = switch off column pivoting
$ pivoting::TRUE $ = switch on column pivoting


Output:

$ Q $ = m-by-m (FULL) or m-by-k (ECONOMY) with k = min(m,n)
$ R $ = m-by-n (FULL) or k-by-n (ECONOMY) with k = min(m,n)
$ p $ = n-element label list
$ P $ = n-by-n permutation matrix

Notes

  • QRMatrix involves modified implementations of the public-domain library TNT, which is available at https://math.nist.gov/tnt/index.html.
  • Q and R are always of the same type of the matrix A.
  • Type can be scalar or complex.
See also
Test-QRMatrix.C
Source files

Definition at line 207 of file QRMatrix.H.

Member Typedef Documentation

◆ cmptType

typedef MatrixType::cmptType cmptType

Definition at line 211 of file QRMatrix.H.

◆ SMatrix

Definition at line 212 of file QRMatrix.H.

◆ RMatrix

Definition at line 213 of file QRMatrix.H.

Member Enumeration Documentation

◆ modes

enum modes : uint8_t

Options for the decomposition mode.

Enumerator
FULL 

compute full-size decomposition

ECONOMY 

compute economy-size decomposition

Definition at line 218 of file QRMatrix.H.

◆ outputs

enum outputs : uint8_t

Options for the output types.

Enumerator
ONLY_Q 

compute only Q

ONLY_R 

compute only R

BOTH_QR 

compute both Q and R

Definition at line 227 of file QRMatrix.H.

◆ pivoting

enum pivoting : bool

Options for the column pivoting.

Enumerator
FALSE 

switch off column pivoting

TRUE 

switch on column pivoting

Definition at line 237 of file QRMatrix.H.

Constructor & Destructor Documentation

◆ QRMatrix() [1/3]

QRMatrix ( )
delete

Construct null.

◆ QRMatrix() [2/3]

QRMatrix ( const modes  mode,
const outputs  output,
const bool  pivoting,
MatrixType &  A 
)
explicit

Construct with a matrix and perform QR decomposition.

Definition at line 376 of file QRMatrix.C.

References A.

◆ QRMatrix() [3/3]

QRMatrix ( const MatrixType &  A,
const modes  mode,
const outputs  output = outputs::BOTH_QR,
const bool  pivoting = false 
)
explicit

Construct with a const matrix and perform QR decomposition.

Definition at line 409 of file QRMatrix.C.

References A.

Member Function Documentation

◆ Q() [1/2]

const MatrixType& Q ( ) const
inlinenoexcept

Return const reference to Q.

Definition at line 381 of file QRMatrix.H.

Referenced by Foam::MatrixTools::pinv().

Here is the caller graph for this function:

◆ Q() [2/2]

MatrixType& Q ( )
inlinenoexcept

Return reference to Q.

Definition at line 389 of file QRMatrix.H.

◆ R() [1/2]

const MatrixType& R ( ) const
inlinenoexcept

Return const reference to R.

Definition at line 397 of file QRMatrix.H.

Referenced by Foam::MatrixTools::pinv().

Here is the caller graph for this function:

◆ R() [2/2]

MatrixType& R ( )
inlinenoexcept

Return reference to R.

Definition at line 405 of file QRMatrix.H.

◆ p()

const labelList& p ( ) const
inlinenoexcept

Return const reference to p.

Definition at line 413 of file QRMatrix.H.

◆ P()

Foam::SquareMatrix< typename MatrixType::cmptType > P ( ) const
inline

Create and return the permutation matrix.

Definition at line 26 of file QRMatrixI.H.

References forAll, and Foam::Zero.

Referenced by Foam::MatrixTools::pinv().

Here is the caller graph for this function:

◆ solve() [1/6]

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

Solve the linear system with the given source and return the solution in the argument x.

Definition at line 442 of file QRMatrix.C.

References x.

Referenced by Foam::MatrixTools::pinv().

Here is the caller graph for this function:

◆ solve() [2/6]

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

Solve the linear system with the given source and return the solution in the argument x.

Definition at line 454 of file QRMatrix.C.

References x.

◆ solve() [3/6]

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

Solve the linear system with the given source and return the solution.

Definition at line 466 of file QRMatrix.C.

◆ solve() [4/6]

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

Solve the linear system with the given source and return the solution.

◆ solve() [5/6]

Foam::RectangularMatrix< typename MatrixType::cmptType > solve ( const RMatrix b)

Solve a row-echelon-form linear system (Ax = b) starting from the bottom by back substitution.

A = R: Non-singular upper-triangular square matrix (m-by-m) b: Source (m-by-p) x: Solution (m-by-p)

Definition at line 497 of file QRMatrix.C.

References Foam::abort(), Foam::constant::atomic::alpha, Foam::constant::physicoChemical::b, Foam::diag(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, k, Matrix< Form, Type >::m(), Foam::mag(), Matrix< Form, Type >::n(), Foam::tab, WarningInFunction, and Foam::Zero.

Here is the call graph for this function:

◆ inv()

Foam::QRMatrix< MatrixType >::SMatrix inv ( ) const

Return the inverse of (Q*R), solving x = (Q*R).inv()*source.

Definition at line 552 of file QRMatrix.C.

References Foam::inv(), Matrix< Form, Type >::m(), and x.

Here is the call graph for this function:

◆ solve() [6/6]

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

Definition at line 482 of file QRMatrix.C.


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