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  Copyright (C) 2024 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 \*---------------------------------------------------------------------------*/
28 
29 #include "scalarMatrices.H"
30 #include "SVD.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
35 (
36  scalarSquareMatrix& matrix,
37  labelList& pivotIndices
38 )
39 {
40  label sign;
41  LUDecompose(matrix, pivotIndices, sign);
42 }
43 
44 
46 (
47  scalarSquareMatrix& matrix,
48  labelList& pivotIndices,
49  label& sign
50 )
51 {
52  const label size = matrix.m();
53 
54  pivotIndices.resize_nocopy(size);
55 
56  List<scalar> vv(size);
57  sign = 1;
58 
59  for (label i = 0; i < size; ++i)
60  {
61  scalar largestCoeff = 0.0;
62  scalar temp;
63  const scalar* __restrict__ matrixi = matrix[i];
64 
65  for (label j = 0; j < size; ++j)
66  {
67  if ((temp = mag(matrixi[j])) > largestCoeff)
68  {
69  largestCoeff = temp;
70  }
71  }
72 
73  if (largestCoeff == 0.0)
74  {
76  << "Singular matrix" << exit(FatalError);
77  }
78 
79  vv[i] = 1.0/largestCoeff;
80  }
81 
82  for (label j = 0; j < size; ++j)
83  {
84  scalar* __restrict__ matrixj = matrix[j];
85 
86  for (label i = 0; i < j; ++i)
87  {
88  scalar* __restrict__ matrixi = matrix[i];
89 
90  scalar sum = matrixi[j];
91  for (label k = 0; k < i; ++k)
92  {
93  sum -= matrixi[k]*matrix(k, j);
94  }
95  matrixi[j] = sum;
96  }
97 
98  label iMax = 0;
99 
100  scalar largestCoeff = 0.0;
101  for (label i = j; i < size; ++i)
102  {
103  scalar* __restrict__ matrixi = matrix[i];
104  scalar sum = matrixi[j];
105 
106  for (label k = 0; k < j; ++k)
107  {
108  sum -= matrixi[k]*matrix(k, j);
109  }
110 
111  matrixi[j] = sum;
112 
113  scalar temp;
114  if ((temp = vv[i]*mag(sum)) >= largestCoeff)
115  {
116  largestCoeff = temp;
117  iMax = i;
118  }
119  }
120 
121  pivotIndices[j] = iMax;
122 
123  if (j != iMax)
124  {
125  scalar* __restrict__ matrixiMax = matrix[iMax];
126 
127  for (label k = 0; k < size; ++k)
128  {
129  std::swap(matrixj[k], matrixiMax[k]);
130  }
131 
132  sign *= -1;
133  vv[iMax] = vv[j];
134  }
135 
136  if (matrixj[j] == 0.0)
137  {
138  matrixj[j] = SMALL;
139  }
140 
141  if (j != size-1)
142  {
143  scalar rDiag = 1.0/matrixj[j];
144 
145  for (label i = j + 1; i < size; ++i)
146  {
147  matrix(i, j) *= rDiag;
148  }
149  }
150  }
151 }
152 
153 
155 {
156  // Store result in upper triangular part of matrix
157  const label size = matrix.m();
158 
159  // Set upper triangular parts to zero.
160  for (label j = 0; j < size; ++j)
161  {
162  for (label k = j + 1; k < size; ++k)
163  {
164  matrix(j, k) = 0.0;
165  }
166  }
167 
168  for (label j = 0; j < size; ++j)
169  {
170  scalar d = 0.0;
171 
172  for (label k = 0; k < j; ++k)
173  {
174  scalar s = 0.0;
175 
176  for (label i = 0; i < k; ++i)
177  {
178  s += matrix(i, k)*matrix(i, j);
179  }
180 
181  s = (matrix(j, k) - s)/matrix(k, k);
182 
183  matrix(k, j) = s;
184  matrix(j, k) = s;
185 
186  d += sqr(s);
187  }
188 
189  d = matrix(j, j) - d;
190 
191  if (d < 0.0)
192  {
194  << "Matrix is not symmetric positive-definite. Unable to "
195  << "decompose."
196  << abort(FatalError);
197  }
198 
199  matrix(j, j) = sqrt(d);
200  }
201 }
202 
203 
204 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
205 
206 void Foam::multiply
207 (
208  scalarRectangularMatrix& ans, // value changed in return
209  const scalarRectangularMatrix& A,
210  const scalarRectangularMatrix& B,
212 )
213 {
214  if (A.n() != B.m())
215  {
217  << "A and B must have identical inner dimensions but A.n = "
218  << A.n() << " and B.m = " << B.m()
219  << abort(FatalError);
220  }
221 
222  if (B.n() != C.m())
223  {
225  << "B and C must have identical inner dimensions but B.n = "
226  << B.n() << " and C.m = " << C.m()
227  << abort(FatalError);
228  }
229 
230  ans = scalarRectangularMatrix(A.m(), C.n(), Foam::zero{});
231 
232  for (label i = 0; i < A.m(); ++i)
233  {
234  for (label g = 0; g < C.n(); ++g)
235  {
236  for (label l = 0; l < C.m(); ++l)
237  {
238  scalar ab = 0;
239  for (label j = 0; j < A.n(); ++j)
240  {
241  ab += A(i, j)*B(j, l);
242  }
243  ans(i, g) += C(l, g) * ab;
244  }
245  }
246  }
247 }
248 
249 
250 void Foam::multiply
251 (
252  scalarRectangularMatrix& ans, // value changed in return
253  const scalarRectangularMatrix& A,
254  const DiagonalMatrix<scalar>& B,
256 )
257 {
258  if (A.n() != B.size())
259  {
261  << "A and B must have identical inner dimensions but A.n = "
262  << A.n() << " and B.m = " << B.size()
263  << abort(FatalError);
264  }
265 
266  if (B.size() != C.m())
267  {
269  << "B and C must have identical inner dimensions but B.n = "
270  << B.size() << " and C.m = " << C.m()
271  << abort(FatalError);
272  }
273 
274  ans = scalarRectangularMatrix(A.m(), C.n(), Zero);
275 
276  for (label i = 0; i < A.m(); ++i)
277  {
278  for (label g = 0; g < C.n(); ++g)
279  {
280  for (label l = 0; l < C.m(); ++l)
281  {
282  ans(i, g) += C(l, g) * A(i, l)*B[l];
283  }
284  }
285  }
286 }
287 
288 
289 void Foam::multiply
290 (
291  scalarSquareMatrix& ans, // value changed in return
292  const scalarSquareMatrix& A,
293  const DiagonalMatrix<scalar>& B,
294  const scalarSquareMatrix& C
295 )
296 {
297  if (A.m() != B.size())
298  {
300  << "A and B must have identical dimensions but A.m = "
301  << A.m() << " and B.m = " << B.size()
302  << abort(FatalError);
303  }
304 
305  if (B.size() != C.m())
306  {
308  << "B and C must have identical dimensions but B.m = "
309  << B.size() << " and C.m = " << C.m()
310  << abort(FatalError);
311  }
312 
313  const label size = A.m();
314 
315  ans = scalarSquareMatrix(size, Zero);
316 
317  for (label i = 0; i < size; ++i)
318  {
319  for (label g = 0; g < size; ++g)
320  {
321  for (label l = 0; l < size; ++l)
322  {
323  ans(i, g) += C(l, g)*A(i, l)*B[l];
324  }
325  }
326  }
327 }
328 
329 
330 //- Pseudo-inverse algorithm for scalar matrices, using Moore-Penrose inverse
332 (
333  const scalarRectangularMatrix& A,
334  scalar minCondition
335 )
336 {
337  SVD svd(A, minCondition);
338  return svd.VSinvUt();
339 }
340 
341 
342 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
scalarRectangularMatrix SVDinv(const scalarRectangularMatrix &A, scalar minCondition=0)
Return the inverse of matrix A using SVD.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const uniformDimensionedVectorField & g
SymmetricSquareMatrix< scalar > scalarSymmetricSquareMatrix
volScalarField & C
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
List< label > labelList
A List of labels.
Definition: List.H:62
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
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))
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
SquareMatrix< scalar > scalarSquareMatrix
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
RectangularMatrix< scalar > scalarRectangularMatrix
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127