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 | Copyright (C) 2011-2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "scalarMatrices.H"
27 #include "SVD.H"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
32 (
33  scalarSquareMatrix& matrix,
34  labelList& pivotIndices
35 )
36 {
37  label sign;
38  LUDecompose(matrix, pivotIndices, sign);
39 }
40 
41 
43 (
44  scalarSquareMatrix& matrix,
45  labelList& pivotIndices,
46  label& sign
47 )
48 {
49  label n = matrix.n();
50  scalar vv[n];
51  sign = 1;
52 
53  for (label i=0; i<n; i++)
54  {
55  scalar largestCoeff = 0.0;
56  scalar temp;
57  const scalar* __restrict__ matrixi = matrix[i];
58 
59  for (label j=0; j<n; j++)
60  {
61  if ((temp = mag(matrixi[j])) > largestCoeff)
62  {
63  largestCoeff = temp;
64  }
65  }
66 
67  if (largestCoeff == 0.0)
68  {
70  (
71  "LUdecompose"
72  "(scalarSquareMatrix& matrix, labelList& rowIndices)"
73  ) << "Singular matrix" << exit(FatalError);
74  }
75 
76  vv[i] = 1.0/largestCoeff;
77  }
78 
79  for (label j=0; j<n; j++)
80  {
81  scalar* __restrict__ matrixj = matrix[j];
82 
83  for (label i=0; i<j; i++)
84  {
85  scalar* __restrict__ matrixi = matrix[i];
86 
87  scalar sum = matrixi[j];
88  for (label k=0; k<i; k++)
89  {
90  sum -= matrixi[k]*matrix[k][j];
91  }
92  matrixi[j] = sum;
93  }
94 
95  label iMax = 0;
96 
97  scalar largestCoeff = 0.0;
98  for (label i=j; i<n; i++)
99  {
100  scalar* __restrict__ matrixi = matrix[i];
101  scalar sum = matrixi[j];
102 
103  for (label k=0; k<j; k++)
104  {
105  sum -= matrixi[k]*matrix[k][j];
106  }
107 
108  matrixi[j] = sum;
109 
110  scalar temp;
111  if ((temp = vv[i]*mag(sum)) >= largestCoeff)
112  {
113  largestCoeff = temp;
114  iMax = i;
115  }
116  }
117 
118  pivotIndices[j] = iMax;
119 
120  if (j != iMax)
121  {
122  scalar* __restrict__ matrixiMax = matrix[iMax];
123 
124  for (label k=0; k<n; k++)
125  {
126  Swap(matrixj[k], matrixiMax[k]);
127  }
128 
129  sign *= -1;
130  vv[iMax] = vv[j];
131  }
132 
133  if (matrixj[j] == 0.0)
134  {
135  matrixj[j] = SMALL;
136  }
137 
138  if (j != n-1)
139  {
140  scalar rDiag = 1.0/matrixj[j];
141 
142  for (label i=j+1; i<n; i++)
143  {
144  matrix[i][j] *= rDiag;
145  }
146  }
147  }
148 }
149 
150 
152 {
153  // Store result in upper triangular part of matrix
154  label size = matrix.n();
155 
156  // Set upper triangular parts to zero.
157  for (label j = 0; j < size; j++)
158  {
159  for (label k = j + 1; k < size; k++)
160  {
161  matrix[j][k] = 0.0;
162  }
163  }
164 
165  for (label j = 0; j < size; j++)
166  {
167  scalar d = 0.0;
168 
169  for (label k = 0; k < j; k++)
170  {
171  scalar s = 0.0;
172 
173  for (label i = 0; i < k; i++)
174  {
175  s += matrix[i][k]*matrix[i][j];
176  }
177 
178  s = (matrix[j][k] - s)/matrix[k][k];
179 
180  matrix[k][j] = s;
181  matrix[j][k] = s;
182 
183  d += sqr(s);
184  }
185 
186  d = matrix[j][j] - d;
187 
188  if (d < 0.0)
189  {
190  FatalErrorIn("Foam::LUDecompose(scalarSymmetricSquareMatrix&)")
191  << "Matrix is not symmetric positive-definite. Unable to "
192  << "decompose."
193  << abort(FatalError);
194  }
195 
196  matrix[j][j] = sqrt(d);
197  }
198 }
199 
200 
201 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
202 
203 void Foam::multiply
204 (
205  scalarRectangularMatrix& ans, // value changed in return
206  const scalarRectangularMatrix& A,
207  const scalarRectangularMatrix& B,
208  const scalarRectangularMatrix& C
209 )
210 {
211  if (A.m() != B.n())
212  {
214  (
215  "multiply("
216  "const scalarRectangularMatrix& A, "
217  "const scalarRectangularMatrix& B, "
218  "const scalarRectangularMatrix& C, "
219  "scalarRectangularMatrix& answer)"
220  ) << "A and B must have identical inner dimensions but A.m = "
221  << A.m() << " and B.n = " << B.n()
222  << abort(FatalError);
223  }
224 
225  if (B.m() != C.n())
226  {
228  (
229  "multiply("
230  "const scalarRectangularMatrix& A, "
231  "const scalarRectangularMatrix& B, "
232  "const scalarRectangularMatrix& C, "
233  "scalarRectangularMatrix& answer)"
234  ) << "B and C must have identical inner dimensions but B.m = "
235  << B.m() << " and C.n = " << C.n()
236  << abort(FatalError);
237  }
238 
239  ans = scalarRectangularMatrix(A.n(), C.m(), scalar(0));
240 
241  for (label i = 0; i < A.n(); i++)
242  {
243  for (label g = 0; g < C.m(); g++)
244  {
245  for (label l = 0; l < C.n(); l++)
246  {
247  scalar ab = 0;
248  for (label j = 0; j < A.m(); j++)
249  {
250  ab += A[i][j]*B[j][l];
251  }
252  ans[i][g] += C[l][g] * ab;
253  }
254  }
255  }
256 }
257 
258 
259 void Foam::multiply
260 (
261  scalarRectangularMatrix& ans, // value changed in return
262  const scalarRectangularMatrix& A,
263  const DiagonalMatrix<scalar>& B,
264  const scalarRectangularMatrix& C
265 )
266 {
267  if (A.m() != B.size())
268  {
270  (
271  "multiply("
272  "const scalarRectangularMatrix& A, "
273  "const DiagonalMatrix<scalar>& B, "
274  "const scalarRectangularMatrix& C, "
275  "scalarRectangularMatrix& answer)"
276  ) << "A and B must have identical inner dimensions but A.m = "
277  << A.m() << " and B.n = " << B.size()
278  << abort(FatalError);
279  }
280 
281  if (B.size() != C.n())
282  {
284  (
285  "multiply("
286  "const scalarRectangularMatrix& A, "
287  "const DiagonalMatrix<scalar>& B, "
288  "const scalarRectangularMatrix& C, "
289  "scalarRectangularMatrix& answer)"
290  ) << "B and C must have identical inner dimensions but B.m = "
291  << B.size() << " and C.n = " << C.n()
292  << abort(FatalError);
293  }
294 
295  ans = scalarRectangularMatrix(A.n(), C.m(), scalar(0));
296 
297  for (label i = 0; i < A.n(); i++)
298  {
299  for (label g = 0; g < C.m(); g++)
300  {
301  for (label l = 0; l < C.n(); l++)
302  {
303  ans[i][g] += C[l][g] * A[i][l]*B[l];
304  }
305  }
306  }
307 }
308 
309 
311 (
312  const scalarRectangularMatrix& A,
313  scalar minCondition
314 )
315 {
316  SVD svd(A, minCondition);
317  return svd.VSinvUt();
318 }
319 
320 
321 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().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 ))
dimensioned< scalar > mag(const dimensioned< Type > &)
RectangularMatrix< scalar > scalarRectangularMatrix
scalarRectangularMatrix SVDinv(const scalarRectangularMatrix &A, scalar minCondition=0)
Return the inverse of matrix A using SVD.
label n() const
Return the number of rows.
Definition: MatrixI.H:56
A templated 2D square symmetric matrix of objects of <T>, where the n x n matrix dimension is known a...
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar sign(const dimensionedScalar &ds)
label n
Singular value decomposition of a rectangular matrix.
Definition: SVD.H:52
label k
Boltzmann constant.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void Swap(T &a, T &b)
Definition: Swap.H:43
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
const dimensionedVector & g
error FatalError
const scalarRectangularMatrix & VSinvUt() const
Return VSinvUt (the pseudo inverse)
Definition: SVDI.H:52
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
label m() const
Return the number of columns.
Definition: MatrixI.H:63