SVD.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) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::SVD
28
29Description
30 Singular value decomposition of a rectangular matrix.
31
32SourceFiles
33 SVDI.H
34 SVD.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef SVD_H
39#define SVD_H
40
41#include "scalarMatrices.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50/*---------------------------------------------------------------------------*\
51 Class SVD Declaration
52\*---------------------------------------------------------------------------*/
54class SVD
55{
56 // Private data
57
58 //- Rectangular matrix with the same dimensions as the input
60
61 //- Square matrix V
63
64 //- The singular values
66
67 //- Convergence flag
68 bool converged_;
69
70 //- The number of zero singular values
71 label nZeros_;
72
73
74 // Private Member Functions
75
76 //- No copy construct
77 SVD(const SVD&) = delete;
78
79 //- No copy assignment
80 void operator=(const SVD&) = delete;
81
82 template<class T>
83 inline const T sign(const T& a, const T& b);
84
85
86public:
87
88 // Constructors
89
90 //- Construct from a rectangular Matrix
91 SVD(const scalarRectangularMatrix& A, const scalar minCondition = 0);
92
93
94 // Access functions
95
96 //- Return U
97 inline const scalarRectangularMatrix& U() const;
98
99 //- Return the square matrix V
100 inline const scalarRectangularMatrix& V() const;
101
102 //- Return the singular values
103 inline const scalarDiagonalMatrix& S() const;
104
105 //- Return the minimum non-zero singular value
106 inline bool converged() const;
107
108 //- Return the number of zero singular values
109 inline label nZeros() const;
110
111 //- Return the minimum non-zero singular value
112 inline scalar minNonZeroS() const;
113
114 //- Return the matrix product V S^(-1) U^T (the pseudo inverse)
116};
117
118
119// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120
121} // End namespace Foam
122
123// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124
125#include "SVDI.H"
126
127// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128
129#endif
130
131// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
A templated (N x N) diagonal matrix of objects of <Type>, effectively containing N elements,...
Singular value decomposition of a rectangular matrix.
Definition: SVD.H:54
const scalarRectangularMatrix & V() const
Return the square matrix V.
Definition: SVDI.H:46
const scalarRectangularMatrix & U() const
Return U.
Definition: SVDI.H:40
scalarRectangularMatrix VSinvUt() const
Return the matrix product V S^(-1) U^T (the pseudo inverse)
Definition: SVD.C:387
scalar minNonZeroS() const
Return the minimum non-zero singular value.
Definition: SVDI.H:70
label nZeros() const
Return the number of zero singular values.
Definition: SVDI.H:64
bool converged() const
Return the minimum non-zero singular value.
Definition: SVDI.H:58
const scalarDiagonalMatrix & S() const
Return the singular values.
Definition: SVDI.H:52
const volScalarField & T
Namespace for OpenFOAM.
volScalarField & b
Definition: createFields.H:27