Loading...
Searching...
No Matches
QRMatrix.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016 OpenFOAM Foundation
9 Copyright (C) 2019-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::QRMatrix
29
30Description
31 \c QRMatrix computes QR decomposition of a given
32 scalar/complex matrix \c A into the following:
33
34 \verbatim
35 A = Q R
36 \endverbatim
37
38 or in case of QR decomposition with column pivoting:
39
40 \verbatim
41 A P = Q R
42 \endverbatim
43
44 where
45 \vartable
46 Q | Unitary/orthogonal matrix
47 R | Upper triangular matrix
48 P | Permutation matrix
49 \endvartable
50
51 References:
52 \verbatim
53 TNT implementation:
54 Pozo, R. (1997).
55 Template Numerical Toolkit for linear algebra:
56 High performance programming with C++
57 and the Standard Template Library.
58 The International Journal of Supercomputer Applications
59 and High Performance Computing, 11(3), 251-263.
60 DOI:10.1177/109434209701100307
61
62 QR decomposition with column pivoting (tag:QSB):
63 Quintana-Ortí, G., Sun, X., & Bischof, C. H. (1998).
64 A BLAS-3 version of the QR factorization with column pivoting.
65 SIAM Journal on Scientific Computing, 19(5), 1486-1494.
66 DOI:10.1137/S1064827595296732
67
68 Moore-Penrose inverse algorithm (tags:KP; KPP):
69 Katsikis, V. N., & Pappas, D. (2008).
70 Fast computing of the Moore-Penrose inverse matrix.
71 Electronic Journal of Linear Algebra, 17(1), 637-650.
72 DOI:10.13001/1081-3810.1287
73
74 Katsikis, V. N., Pappas, D., & Petralias, A. (2011).
75 An improved method for the computation of
76 the Moore–Penrose inverse matrix.
77 Applied Mathematics and Computation, 217(23), 9828-9834.
78 DOI:10.1016/j.amc.2011.04.080
79
80 Tolerance for the Moore-Penrose inverse algorithm (tag:TA):
81 Toutounian, F., & Ataei, A. (2009).
82 A new method for computing Moore–Penrose inverse matrices.
83 Journal of Computational and applied Mathematics, 228(1), 412-417.
84 DOI:10.1016/j.cam.2008.10.008
85 \endverbatim
86
87Usage
88
89 \heading Input:
90 \vartable
91 A | RectangularMatrix<Type> or SquareMatrix<Type>
92 \endvartable
93
94 \heading Options for the decomposition mode:
95 \vartable
96 modes::FULL | compute full-size decomposition
97 modes::ECONOMY | compute economy-size decomposition
98 \endvartable
99
100 \heading Options for the output types:
101 \vartable
102 outputs::ONLY\_Q | compute only Q
103 outputs::ONLY\_R | compute only R
104 outputs::BOTH\_QR | compute both Q and R
105 \endvartable
106
107 \heading Options for the column pivoting:
108 \vartable
109 pivoting::FALSE | switch off column pivoting
110 pivoting::TRUE | switch on column pivoting
111 \endvartable
112
113 \heading Output:
114 \vartable
115 Q | m-by-m (FULL) or m-by-k (ECONOMY) with k = min(m,n)
116 R | m-by-n (FULL) or k-by-n (ECONOMY) with k = min(m,n)
117 p | n-element label list
118 P | n-by-n permutation matrix
119 \endvartable
120
121Notes
122 - \c QRMatrix involves modified implementations of the public-domain
123 library \c TNT, which is available at https://math.nist.gov/tnt/index.html.
124 - \c Q and \c R are always of the same type of the matrix \c A.
125 - \c Type can be \c scalar or \c complex.
126
127See also
128 Test-QRMatrix.C
129
130SourceFiles
131 QRMatrix.C
132 QRMatrixI.H
133
134\*---------------------------------------------------------------------------*/
135
136#ifndef Foam_QRMatrix_H
137#define Foam_QRMatrix_H
138
139#include "RectangularMatrix.H"
140#include "SquareMatrix.H"
141
142// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143
144namespace Foam
145{
146
147/*---------------------------------------------------------------------------*\
148 Class QRMatrix Declaration
149\*---------------------------------------------------------------------------*/
150
151template<class MatrixType>
152class QRMatrix
153{
154public:
155
156 typedef typename MatrixType::cmptType cmptType;
157 typedef typename MatrixType::size_type size_type;
158 typedef typename MatrixType::value_type value_type;
159
160 typedef SquareMatrix<cmptType> SMatrix;
161 typedef RectangularMatrix<cmptType> RMatrix;
162
163 //- Options for the decomposition mode
164 enum modes : uint8_t
165 {
166 FULL = 1,
167 ECONOMY = 2,
168 };
169
170 //- Options for the output types
171 enum outputs : uint8_t
172 {
173 ONLY_Q = 1,
174 ONLY_R = 2,
175 BOTH_QR = 3
176 };
177
178 //- Options for the column pivoting
179 enum pivoting : bool
180 {
181 FALSE = false,
182 TRUE = true
183 };
184
185
186private:
187
188 // Private Data
189
190 //- Decomposition mode
191 const modes mode_;
192
193 //- Output type
194 const outputs output_;
195
196 //- Output shape factor based on the decomposition mode
197 const label sz_;
198
199 //- Unitary/orthogonal matrix
200 MatrixType Q_;
201
202 //- Upper triangular matrix
203 MatrixType R_;
204
205 //- Permutation list (if column-pivoting is on)
206 labelList p_;
208
209 // Private Member Functions
210
211 // Evaluation
213 //- Calculate output shape factor based on the decomposition mode
214 label calcShapeFactor(const MatrixType& A) const;
216 //- Compute QR decomposition
217 void decompose(MatrixType& AT);
218
219 //- Compute QR decomposition with column pivoting
220 void decompose(MatrixType& AT, const bool pivot);
222 //- Calculate Q
223 void calcQ(const MatrixType& AT);
225 //- Calculate R
226 void calcR(const MatrixType& AT, const List<cmptType>& diag);
227
228
229 // Linear system solution
231 //- Solve the linear system with the Field argument x
232 //- initialized to the appropriate transformed source
233 //- (e.g. Q.T()*source) and return the solution in x
234 template<template<typename> class ListContainer>
235 void solvex(ListContainer<cmptType>& x) const;
236
237 //- Solve the linear system with the given source
238 //- and return the solution in x
239 template<template<typename> class ListContainer>
240 void solveImpl
241 (
243 const ListContainer<cmptType>& source
244 ) const;
245
246
247public:
248
249 // Generated Methods
250
251 //- No default construct
252 QRMatrix() = delete;
253
254 //- No copy construct
255 QRMatrix(const QRMatrix&) = delete;
256
257 //- No copy assignment
258 QRMatrix& operator=(const QRMatrix&) = delete;
259
260
261 // Constructors
262
263 //- Construct with a matrix and perform QR decomposition
264 explicit QRMatrix
265 (
266 const modes mode,
267 const outputs output,
268 const bool pivoting,
269 MatrixType& A
270 );
271
272 //- Construct with a const matrix and perform QR decomposition
273 explicit QRMatrix
274 (
275 const MatrixType& A,
276 const modes mode,
277 const outputs output = outputs::BOTH_QR,
278 const bool pivoting = false
279 );
280
281
282 // Member Functions
283
284 // Access
285
286 //- Return const reference to Q
287 const MatrixType& Q() const noexcept
288 {
289 return Q_;
290 }
291
292 //- Return reference to Q
293 MatrixType& Q() noexcept
294 {
295 return Q_;
296 }
297
298 //- Return const reference to R
299 const MatrixType& R() const noexcept
300 {
301 return R_;
302 }
303
304 //- Return reference to R
305 MatrixType& R() noexcept
306 {
307 return R_;
308 }
309
310 //- Return const reference to p
311 const labelList& p() const noexcept
312 {
313 return p_;
314 }
315
316 //- Create and return the permutation matrix
317 inline SMatrix P() const;
318
319
320 // Linear system solution
321
322 //- Solve the linear system with the given source
323 //- and return the solution in the argument x
324 void solve
325 (
326 List<cmptType>& x,
327 const UList<cmptType>& source
328 ) const;
329
330 //- Solve the linear system with the given source
331 //- and return the solution in the argument x
332 template<class Addr>
333 void solve
334 (
335 List<cmptType>& x,
336 const IndirectListBase<cmptType, Addr>& source
337 ) const;
338
339 //- Solve the linear system with the given source
340 //- and return the solution
342 (
343 const UList<cmptType>& source
344 ) const;
345
346 //- Solve the linear system with the given source
347 //- and return the solution
348 template<class Addr>
350 (
352 ) const;
353
354 //- Solve a row-echelon-form linear system (Ax = b)
355 //- starting from the bottom by back substitution
356 // A = R: Non-singular upper-triangular square matrix (m-by-m)
357 // b: Source (m-by-p)
358 // x: Solution (m-by-p)
360 (
361 const RMatrix& b
362 );
363
364 //- Return the inverse of (Q*R), solving x = (Q*R).inv()*source
365 SMatrix inv() const;
366};
367
368
369// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
370
371namespace MatrixTools
372{
373
374// * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * * //
375
376//- Moore-Penrose inverse of singular/non-singular square/rectangular
377//- scalar/complex matrices (KPP:p. 9834; KP:p. 648)
378// The tolerance to ensure the R1 matrix full-rank is set to 1e-5
379// by (TA; mentioned in (KPP:p. 9832)) in contrast to 1e-13 (KPP:p. 9834).
380template<class MatrixType>
381MatrixType pinv
382(
383 const MatrixType& A,
384 scalar tol = 1e-5
385);
387// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388
389} // End namespace MatrixTools
390} // End namespace Foam
391
392// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393
394#include "QRMatrixI.H"
395
396// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
397
398#ifdef NoRepository
399 #include "QRMatrix.C"
400#endif
401
402// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403
404#endif
405
406// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Base for lists with indirect addressing, templated on the list contents type and the addressing type....
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
QRMatrix computes QR decomposition of a given scalar/complex matrix A into the following:
Definition QRMatrix.H:208
tmp< Field< cmptType > > solve(const IndirectListBase< cmptType, Addr > &source) const
Solve the linear system with the given source and return the solution.
MatrixType & Q() noexcept
Return reference to Q.
Definition QRMatrix.H:394
SMatrix inv() const
Return the inverse of (Q*R), solving x = (Q*R).inv()*source.
Definition QRMatrix.C:552
QRMatrix(const QRMatrix &)=delete
No copy construct.
MatrixType::size_type size_type
Definition QRMatrix.H:212
MatrixType::value_type value_type
Definition QRMatrix.H:213
QRMatrix()=delete
No default construct.
pivoting
Options for the column pivoting.
Definition QRMatrix.H:241
@ FALSE
switch off column pivoting
Definition QRMatrix.H:242
@ TRUE
switch on column pivoting
Definition QRMatrix.H:243
SMatrix P() const
Create and return the permutation matrix.
Definition QRMatrixI.H:26
outputs
Options for the output types.
Definition QRMatrix.H:231
@ ONLY_Q
compute only Q
Definition QRMatrix.H:232
@ ONLY_R
compute only R
Definition QRMatrix.H:233
@ BOTH_QR
compute both Q and R
Definition QRMatrix.H:234
modes
Options for the decomposition mode.
Definition QRMatrix.H:222
@ ECONOMY
compute economy-size decomposition
Definition QRMatrix.H:224
@ FULL
compute full-size decomposition
Definition QRMatrix.H:223
SquareMatrix< cmptType > SMatrix
Definition QRMatrix.H:215
RectangularMatrix< cmptType > RMatrix
Definition QRMatrix.H:216
MatrixType::cmptType cmptType
Definition QRMatrix.H:211
const MatrixType & Q() const noexcept
Return const reference to Q.
Definition QRMatrix.H:386
MatrixType & R() noexcept
Return reference to R.
Definition QRMatrix.H:410
const labelList & p() const noexcept
Return const reference to p.
Definition QRMatrix.H:418
const MatrixType & R() const noexcept
Return const reference to R.
Definition QRMatrix.H:402
QRMatrix & operator=(const QRMatrix &)=delete
No copy assignment.
A templated (M x N) rectangular matrix of objects of <Type>, containing M*N elements,...
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
A class for managing temporary objects.
Definition tmp.H:75
Collection of functions for matrix-related verifications.
MatrixType pinv(const MatrixType &A, scalar tol=1e-5)
Moore-Penrose inverse of singular/non-singular square/rectangular scalar/complex matrices (KPP:p....
Definition QRMatrix.C:577
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
const direction noexcept
Definition scalarImpl.H:265
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
CEqn solve()
volScalarField & b
volScalarField & e