scalarMatrices.C
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
26\*---------------------------------------------------------------------------*/
27
28#include "scalarMatrices.H"
29#include "SVD.H"
30
31// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32
34(
35 scalarSquareMatrix& matrix,
36 labelList& pivotIndices
37)
38{
39 label sign;
40 LUDecompose(matrix, pivotIndices, sign);
41}
42
43
45(
46 scalarSquareMatrix& matrix,
47 labelList& pivotIndices,
48 label& sign
49)
50{
51 label m = matrix.m();
52 scalar vv[m];
53 sign = 1;
54
55 for (label i = 0; i < m; ++i)
56 {
57 scalar largestCoeff = 0.0;
58 scalar temp;
59 const scalar* __restrict__ matrixi = matrix[i];
60
61 for (label j = 0; j < m; ++j)
62 {
63 if ((temp = mag(matrixi[j])) > largestCoeff)
64 {
65 largestCoeff = temp;
66 }
67 }
68
69 if (largestCoeff == 0.0)
70 {
72 << "Singular matrix" << exit(FatalError);
73 }
74
75 vv[i] = 1.0/largestCoeff;
76 }
77
78 for (label j = 0; j < m; ++j)
79 {
80 scalar* __restrict__ matrixj = matrix[j];
81
82 for (label i = 0; i < j; ++i)
83 {
84 scalar* __restrict__ matrixi = matrix[i];
85
86 scalar sum = matrixi[j];
87 for (label k = 0; k < i; ++k)
88 {
89 sum -= matrixi[k]*matrix(k, j);
90 }
91 matrixi[j] = sum;
92 }
93
94 label iMax = 0;
95
96 scalar largestCoeff = 0.0;
97 for (label i = j; i < m; ++i)
98 {
99 scalar* __restrict__ matrixi = matrix[i];
100 scalar sum = matrixi[j];
101
102 for (label k = 0; k < j; ++k)
103 {
104 sum -= matrixi[k]*matrix(k, j);
105 }
106
107 matrixi[j] = sum;
108
109 scalar temp;
110 if ((temp = vv[i]*mag(sum)) >= largestCoeff)
111 {
112 largestCoeff = temp;
113 iMax = i;
114 }
115 }
116
117 pivotIndices[j] = iMax;
118
119 if (j != iMax)
120 {
121 scalar* __restrict__ matrixiMax = matrix[iMax];
122
123 for (label k = 0; k < m; ++k)
124 {
125 Swap(matrixj[k], matrixiMax[k]);
126 }
127
128 sign *= -1;
129 vv[iMax] = vv[j];
130 }
131
132 if (matrixj[j] == 0.0)
133 {
134 matrixj[j] = SMALL;
135 }
136
137 if (j != m-1)
138 {
139 scalar rDiag = 1.0/matrixj[j];
140
141 for (label i = j + 1; i < m; ++i)
142 {
143 matrix(i, j) *= rDiag;
144 }
145 }
146 }
147}
148
149
151{
152 // Store result in upper triangular part of matrix
153 label size = matrix.m();
154
155 // Set upper triangular parts to zero.
156 for (label j = 0; j < size; ++j)
157 {
158 for (label k = j + 1; k < size; ++k)
159 {
160 matrix(j, k) = 0.0;
161 }
162 }
163
164 for (label j = 0; j < size; ++j)
165 {
166 scalar d = 0.0;
167
168 for (label k = 0; k < j; ++k)
169 {
170 scalar s = 0.0;
171
172 for (label i = 0; i < k; ++i)
173 {
174 s += matrix(i, k)*matrix(i, j);
175 }
176
177 s = (matrix(j, k) - s)/matrix(k, k);
178
179 matrix(k, j) = s;
180 matrix(j, k) = s;
181
182 d += sqr(s);
183 }
184
185 d = matrix(j, j) - d;
186
187 if (d < 0.0)
188 {
190 << "Matrix is not symmetric positive-definite. Unable to "
191 << "decompose."
192 << abort(FatalError);
193 }
194
195 matrix(j, j) = sqrt(d);
196 }
197}
198
199
200// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
201
203(
204 scalarRectangularMatrix& ans, // value changed in return
208)
209{
210 if (A.n() != B.m())
211 {
213 << "A and B must have identical inner dimensions but A.n = "
214 << A.n() << " and B.m = " << B.m()
215 << abort(FatalError);
216 }
217
218 if (B.n() != C.m())
219 {
221 << "B and C must have identical inner dimensions but B.n = "
222 << B.n() << " and C.m = " << C.m()
223 << abort(FatalError);
224 }
225
226 ans = scalarRectangularMatrix(A.m(), C.n(), Zero);
227
228 for (label i = 0; i < A.m(); ++i)
229 {
230 for (label g = 0; g < C.n(); ++g)
231 {
232 for (label l = 0; l < C.m(); ++l)
233 {
234 scalar ab = 0;
235 for (label j = 0; j < A.n(); ++j)
236 {
237 ab += A(i, j)*B(j, l);
238 }
239 ans(i, g) += C(l, g) * ab;
240 }
241 }
242 }
243}
244
245
247(
248 scalarRectangularMatrix& ans, // value changed in return
252)
253{
254 if (A.n() != B.size())
255 {
257 << "A and B must have identical inner dimensions but A.n = "
258 << A.n() << " and B.m = " << B.size()
259 << abort(FatalError);
260 }
261
262 if (B.size() != C.m())
263 {
265 << "B and C must have identical inner dimensions but B.n = "
266 << B.size() << " and C.m = " << C.m()
267 << abort(FatalError);
268 }
269
270 ans = scalarRectangularMatrix(A.m(), C.n(), Zero);
271
272 for (label i = 0; i < A.m(); ++i)
273 {
274 for (label g = 0; g < C.n(); ++g)
275 {
276 for (label l = 0; l < C.m(); ++l)
277 {
278 ans(i, g) += C(l, g) * A(i, l)*B[l];
279 }
280 }
281 }
282}
283
284
286(
287 scalarSquareMatrix& ans, // value changed in return
288 const scalarSquareMatrix& A,
290 const scalarSquareMatrix& C
291)
292{
293 if (A.m() != B.size())
294 {
296 << "A and B must have identical dimensions but A.m = "
297 << A.m() << " and B.m = " << B.size()
298 << abort(FatalError);
299 }
300
301 if (B.size() != C.m())
302 {
304 << "B and C must have identical dimensions but B.m = "
305 << B.size() << " and C.m = " << C.m()
306 << abort(FatalError);
307 }
308
309 const label size = A.m();
310
311 ans = scalarSquareMatrix(size, Zero);
312
313 for (label i = 0; i < size; ++i)
314 {
315 for (label g = 0; g < size; ++g)
316 {
317 for (label l = 0; l < size; ++l)
318 {
319 ans(i, g) += C(l, g)*A(i, l)*B[l];
320 }
321 }
322 }
323}
324
325
326//- Pseudo-inverse algorithm for scalar matrices, using Moore-Penrose inverse
328(
330 scalar minCondition
331)
332{
333 SVD svd(A, minCondition);
334 return svd.VSinvUt();
335}
336
337
338// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
label k
const uniformDimensionedVectorField & g
Graphite solid properties.
Definition: C.H:53
A templated (N x N) diagonal matrix of objects of <Type>, effectively containing N elements,...
label m() const noexcept
The number of rows.
Definition: MatrixI.H:96
Singular value decomposition of a rectangular matrix.
Definition: SVD.H:54
scalarRectangularMatrix VSinvUt() const
Return the matrix product V S^(-1) U^T (the pseudo inverse)
Definition: SVD.C:387
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
void Swap(DynamicList< T, SizeMinA > &a, DynamicList< T, SizeMinB > &b)
Definition: DynamicList.H:408
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalarRectangularMatrix SVDinv(const scalarRectangularMatrix &A, scalar minCondition=0)
Return the inverse of matrix A using SVD.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130