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 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::QRMatrix
29 
30 Description
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 
87 Usage
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 
121 Notes
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 
127 See also
128  Test-QRMatrix.C
129 
130 SourceFiles
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 
144 namespace Foam
145 {
146 
147 /*---------------------------------------------------------------------------*\
148  Class QRMatrix Declaration
149 \*---------------------------------------------------------------------------*/
150 
151 template<class MatrixType>
152 class QRMatrix
153 {
154 public:
155 
156  typedef typename MatrixType::cmptType cmptType;
157  typedef SquareMatrix<cmptType> SMatrix;
158  typedef RectangularMatrix<cmptType> RMatrix;
159 
160  //- Options for the decomposition mode
161  enum modes : uint8_t
162  {
163  FULL = 1,
164  ECONOMY = 2,
165  };
166 
167  //- Options for the output types
168  enum outputs : uint8_t
169  {
170  ONLY_Q = 1,
171  ONLY_R = 2,
172  BOTH_QR = 3
173  };
174 
175  //- Options for the column pivoting
176  enum pivoting : bool
177  {
178  FALSE = false,
179  TRUE = true
180  };
181 
182 
183 private:
184 
185  // Private Data
186 
187  //- Decomposition mode
188  const modes mode_;
189 
190  //- Output type
191  const outputs output_;
192 
193  //- Output shape factor based on the decomposition mode
194  const label sz_;
195 
196  //- Unitary/orthogonal matrix
197  MatrixType Q_;
198 
199  //- Upper triangular matrix
200  MatrixType R_;
201 
202  //- Permutation list (if column-pivoting is on)
203  labelList p_;
204 
205 
206  // Private Member Functions
208  // Evaluation
209 
210  //- Calculate output shape factor based on the decomposition mode
211  label calcShapeFactor(const MatrixType& A) const;
213  //- Compute QR decomposition
214  void decompose(MatrixType& AT);
215 
216  //- Compute QR decomposition with column pivoting
217  void decompose(MatrixType& AT, const bool pivot);
219  //- Calculate Q
220  void calcQ(const MatrixType& AT);
222  //- Calculate R
223  void calcR(const MatrixType& AT, const List<cmptType>& diag);
224 
225 
226  // Linear system solution
228  //- Solve the linear system with the Field argument x
229  //- initialized to the appropriate transformed source
230  //- (e.g. Q.T()*source) and return the solution in x
231  template<template<typename> class ListContainer>
232  void solvex(ListContainer<cmptType>& x) const;
233 
234  //- Solve the linear system with the given source
235  //- and return the solution in x
236  template<template<typename> class ListContainer>
237  void solveImpl
238  (
240  const ListContainer<cmptType>& source
241  ) const;
242 
243 
244  //- No copy construct
245  QRMatrix(const QRMatrix&) = delete;
246 
247  //- No copy assignment
248  QRMatrix& operator=(const QRMatrix&) = delete;
249 
250 
251 public:
252 
253  // Constructors
254 
255  //- Construct null
256  QRMatrix() = delete;
257 
258  //- Construct with a matrix and perform QR decomposition
259  explicit QRMatrix
260  (
261  const modes mode,
262  const outputs output,
263  const bool pivoting,
264  MatrixType& A
265  );
266 
267  //- Construct with a const matrix and perform QR decomposition
268  explicit QRMatrix
269  (
270  const MatrixType& A,
271  const modes mode,
272  const outputs output = outputs::BOTH_QR,
273  const bool pivoting = false
274  );
275 
276 
277  // Member Functions
278 
279  // Access
280 
281  //- Return const reference to Q
282  const MatrixType& Q() const noexcept
283  {
284  return Q_;
285  }
286 
287  //- Return reference to Q
288  MatrixType& Q() noexcept
289  {
290  return Q_;
291  }
292 
293  //- Return const reference to R
294  const MatrixType& R() const noexcept
295  {
296  return R_;
297  }
298 
299  //- Return reference to R
300  MatrixType& R() noexcept
301  {
302  return R_;
303  }
304 
305  //- Return const reference to p
306  const labelList& p() const noexcept
307  {
308  return p_;
309  }
310 
311  //- Create and return the permutation matrix
312  inline SMatrix P() const;
313 
314 
315  // Linear system solution
316 
317  //- Solve the linear system with the given source
318  //- and return the solution in the argument x
319  void solve
320  (
321  List<cmptType>& x,
322  const UList<cmptType>& source
323  ) const;
324 
325  //- Solve the linear system with the given source
326  //- and return the solution in the argument x
327  template<class Addr>
328  void solve
329  (
330  List<cmptType>& x,
331  const IndirectListBase<cmptType, Addr>& source
332  ) const;
333 
334  //- Solve the linear system with the given source
335  //- and return the solution
336  tmp<Field<cmptType>> solve
337  (
338  const UList<cmptType>& source
339  ) const;
340 
341  //- Solve the linear system with the given source
342  //- and return the solution
343  template<class Addr>
344  tmp<Field<cmptType>> solve
345  (
346  const IndirectListBase<cmptType, Addr>& source
347  ) const;
348 
349  //- Solve a row-echelon-form linear system (Ax = b)
350  //- starting from the bottom by back substitution
351  // A = R: Non-singular upper-triangular square matrix (m-by-m)
352  // b: Source (m-by-p)
353  // x: Solution (m-by-p)
354  RMatrix solve
355  (
356  const RMatrix& b
357  );
358 
359  //- Return the inverse of (Q*R), solving x = (Q*R).inv()*source
360  SMatrix inv() const;
361 };
362 
363 namespace MatrixTools
364 {
365 
366 // * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * * //
367 
368 //- Moore-Penrose inverse of singular/non-singular square/rectangular
369 //- scalar/complex matrices (KPP:p. 9834; KP:p. 648)
370 // The tolerance to ensure the R1 matrix full-rank is set to 1e-5
371 // by (TA; mentioned in (KPP:p. 9832)) in contrast to 1e-13 (KPP:p. 9834).
372 template<class MatrixType>
373 MatrixType pinv
374 (
375  const MatrixType& A,
376  scalar tol = 1e-5
377 );
378 
379 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
380 
381 } // End namespace MatrixTools
382 } // End namespace Foam
383 
384 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
385 
386 #include "QRMatrixI.H"
387 
388 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
390 #ifdef NoRepository
391  #include "QRMatrix.C"
392 #endif
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 #endif
398 // ************************************************************************* //
MatrixType::cmptType cmptType
Definition: QRMatrix.H:211
pivoting
Options for the column pivoting.
Definition: QRMatrix.H:237
outputs
Options for the output types.
Definition: QRMatrix.H:227
RectangularMatrix< cmptType > RMatrix
Definition: QRMatrix.H:213
QRMatrix computes QR decomposition of a given scalar/complex matrix A into the following: ...
Definition: QRMatrix.H:207
SMatrix P() const
Create and return the permutation matrix.
Definition: QRMatrixI.H:26
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: QRMatrix.C:442
SMatrix inv() const
Return the inverse of (Q*R), solving x = (Q*R).inv()*source.
Definition: QRMatrix.C:552
compute economy-size decomposition
Definition: QRMatrix.H:221
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
const MatrixType & R() const noexcept
Return const reference to R.
Definition: QRMatrix.H:397
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
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
compute only Q
Definition: QRMatrix.H:229
const MatrixType & Q() const noexcept
Return const reference to Q.
Definition: QRMatrix.H:381
compute full-size decomposition
Definition: QRMatrix.H:220
const direction noexcept
Definition: Scalar.H:258
const labelList & p() const noexcept
Return const reference to p.
Definition: QRMatrix.H:413
switch off column pivoting
Definition: QRMatrix.H:239
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:44
compute only R
Definition: QRMatrix.H:230
List< label > labelList
A List of labels.
Definition: List.H:62
SquareMatrix< cmptType > SMatrix
Definition: QRMatrix.H:212
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
modes
Options for the decomposition mode.
Definition: QRMatrix.H:218
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: POSIX.C:773
switch on column pivoting
Definition: QRMatrix.H:240
QRMatrix()=delete
Construct null.
compute both Q and R
Definition: QRMatrix.H:231
Namespace for OpenFOAM.