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< cmptType > | SMatrix |
typedef RectangularMatrix< cmptType > | RMatrix |
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 labelList & | p () 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 |
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 |
RMatrix | solve (const RMatrix &b) |
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 |
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
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
.Definition at line 211 of file QRMatrix.H.
typedef MatrixType::cmptType cmptType |
Definition at line 215 of file QRMatrix.H.
typedef SquareMatrix<cmptType> SMatrix |
Definition at line 216 of file QRMatrix.H.
typedef RectangularMatrix<cmptType> RMatrix |
Definition at line 217 of file QRMatrix.H.
enum modes : uint8_t |
Options for the decomposition mode.
Enumerator | |
---|---|
FULL | compute full-size decomposition |
ECONOMY | compute economy-size decomposition |
Definition at line 220 of file QRMatrix.H.
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.
Options for the column pivoting.
Enumerator | |
---|---|
FALSE | switch off column pivoting |
TRUE | switch on column pivoting |
Definition at line 235 of file QRMatrix.H.
|
delete |
Construct null.
Construct with a matrix and perform QR decomposition.
Definition at line 382 of file QRMatrix.C.
References A, and Foam::mode().
|
explicit |
Construct with a const matrix and perform QR decomposition.
Definition at line 415 of file QRMatrix.C.
References A.
|
inlinenoexcept |
Return const reference to Q.
Definition at line 341 of file QRMatrix.H.
Referenced by Foam::MatrixTools::pinv().
|
inlinenoexcept |
Return reference to Q.
Definition at line 347 of file QRMatrix.H.
|
inlinenoexcept |
Return const reference to R.
Definition at line 353 of file QRMatrix.H.
Referenced by Foam::MatrixTools::pinv().
|
inlinenoexcept |
Return reference to R.
Definition at line 359 of file QRMatrix.H.
|
inlinenoexcept |
Return const reference to p.
Definition at line 365 of file QRMatrix.H.
|
inline |
Create and return the permutation matrix.
Definition at line 33 of file QRMatrixI.H.
References forAll.
Referenced by Foam::MatrixTools::pinv().
Solve the linear system with the given source and return the solution in the argument x
Definition at line 448 of file QRMatrix.C.
References x.
Referenced by Foam::MatrixTools::pinv().
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 460 of file QRMatrix.C.
References x.
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 472 of file QRMatrix.C.
tmp< Field< cmptType > > solve | ( | const IndirectListBase< cmptType, Addr > & | source | ) | const |
Solve the linear system with the given source and return the solution
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 503 of file QRMatrix.C.
References Foam::abort(), alpha, 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.
Foam::QRMatrix< MatrixType >::SMatrix inv |
Return the inverse of (Q*R), solving x = (Q*R).inv()*source.
Definition at line 559 of file QRMatrix.C.
References Foam::inv(), Matrix< Form, Type >::m(), and x.
Foam::tmp< Foam::Field< typename MatrixType::cmptType > > solve | ( | const IndirectListBase< cmptType, Addr > & | source | ) | const |
Definition at line 488 of file QRMatrix.C.