All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SVD.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::SVD
26 
27 Description
28  Singular value decomposition of a rectangular matrix.
29 
30 SourceFiles
31  SVDI.H
32  SVD.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef SVD_H
37 #define SVD_H
38 
39 #include "scalarMatrices.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 /*---------------------------------------------------------------------------*\
49  Class SVD Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 class SVD
53 {
54  // Private Data
55 
56  //- Rectangular matrix with the same dimensions as the input
58 
59  //- Square matrix V
61 
62  //- The singular values
64 
65  //- Convergence flag
66  bool converged_;
67 
68  //- The number of zero singular values
69  label nZeros_;
70 
71 
72  // Private Member Functions
73 
74  template<class T>
75  inline const T sign(const T& a, const T& b);
76 
77 
78 public:
79 
80  // Constructors
81 
82  //- Construct from a rectangular Matrix
83  SVD(const scalarRectangularMatrix& A, const scalar minCondition = 0);
84 
85  //- Disallow default bitwise copy construction
86  SVD(const SVD&) = delete;
87 
88 
89  // Access functions
90 
91  //- Return U
92  inline const scalarRectangularMatrix& U() const;
93 
94  //- Return the square matrix V
95  inline const scalarRectangularMatrix& V() const;
96 
97  //- Return the singular values
98  inline const scalarDiagonalMatrix& S() const;
99 
100  //- Return the minimum non-zero singular value
101  inline bool converged() const;
102 
103  //- Return the number of zero singular values
104  inline label nZeros() const;
105 
106  //- Return the minimum non-zero singular value
107  inline scalar minNonZeroS() const;
108 
109  //- Return the matrix product V S^(-1) U^T (the pseudo inverse)
111 
112 
113  // Member Operators
114 
115  //- Disallow default bitwise assignment
116  void operator=(const SVD&) = delete;
117 };
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 } // End namespace Foam
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 #include "SVDI.H"
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
130 #endif
131 
132 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
SVD(const scalarRectangularMatrix &A, const scalar minCondition=0)
Construct from a rectangular Matrix.
Definition: SVD.C:33
label nZeros() const
Return the number of zero singular values.
Definition: SVDI.H:62
const scalarRectangularMatrix & U() const
Return U.
Definition: SVDI.H:38
scalar minNonZeroS() const
Return the minimum non-zero singular value.
Definition: SVDI.H:68
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
const scalarDiagonalMatrix & S() const
Return the singular values.
Definition: SVDI.H:50
bool converged() const
Return the minimum non-zero singular value.
Definition: SVDI.H:56
void operator=(const SVD &)=delete
Disallow default bitwise assignment.
Singular value decomposition of a rectangular matrix.
Definition: SVD.H:51
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalarRectangularMatrix VSinvUt() const
Return the matrix product V S^(-1) U^T (the pseudo inverse)
Definition: SVD.C:385
const scalarRectangularMatrix & V() const
Return the square matrix V.
Definition: SVDI.H:44
Namespace for OpenFOAM.